MaCh3  2.2.3
Reference Guide
ParameterHandlerBase.cpp
Go to the documentation of this file.
2 #include "regex"
3 
4 // ********************************************
5 ParameterHandlerBase::ParameterHandlerBase(std::string name, std::string file, double threshold, int FirstPCA, int LastPCA)
6  : inputFile(file), pca(false) {
7 // ********************************************
8  MACH3LOG_DEBUG("Constructing instance of ParameterHandler");
9  doSpecialStepProposal = false;
10  if (threshold < 0 || threshold >= 1) {
11  MACH3LOG_INFO("NOTE: {} {}", name, file);
12  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");
13  MACH3LOG_INFO("Please specify a number between 0 and 1");
14  MACH3LOG_INFO("You specified: ");
15  MACH3LOG_INFO("Am instead calling the usual non-PCA constructor...");
16  pca = false;
17  }
18  Init(name, file);
19 
20  // Call the innocent helper function
21  if (pca) ConstructPCA(threshold, FirstPCA, LastPCA);
22 }
23 // ********************************************
24 ParameterHandlerBase::ParameterHandlerBase(const std::vector<std::string>& YAMLFile, std::string name, double threshold, int FirstPCA, int LastPCA)
25  : inputFile(YAMLFile[0].c_str()), matrixName(name), pca(true) {
26 // ********************************************
27  MACH3LOG_INFO("Constructing instance of ParameterHandler using");
28  doSpecialStepProposal = false;
29  for(unsigned int i = 0; i < YAMLFile.size(); i++)
30  {
31  MACH3LOG_INFO("{}", YAMLFile[i]);
32  }
33  MACH3LOG_INFO("as an input");
34 
35  if (threshold < 0 || threshold >= 1) {
36  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");
37  MACH3LOG_INFO("Please specify a number between 0 and 1");
38  MACH3LOG_INFO("You specified: ");
39  MACH3LOG_INFO("Am instead calling the usual non-PCA constructor...");
40  pca = false;
41  }
42 
43  Init(YAMLFile);
44  // Call the innocent helper function
45  if (pca) ConstructPCA(threshold, FirstPCA, LastPCA);
46 }
47 
48 // ********************************************
49 //Destructor
51 // ********************************************
52  delete[] randParams;
53  delete[] corr_throw;
54 
55  if (covMatrix != nullptr) delete covMatrix;
56  if (invCovMatrix != nullptr) delete invCovMatrix;
57  if (throwMatrix != nullptr) delete throwMatrix;
58  for(int i = 0; i < _fNumPar; i++) {
59  delete[] throwMatrixCholDecomp[i];
60  }
61  delete[] throwMatrixCholDecomp;
62 }
63 
64 // ********************************************
65 void ParameterHandlerBase::ConstructPCA(const double eigen_threshold, int FirstPCAdpar, int LastPCAdpar) {
66 // ********************************************
67  if(AdaptiveHandler) {
68  MACH3LOG_ERROR("Adaption has been enabled and now trying to enable PCA. Right now both configuration don't work with each other");
69  throw MaCh3Exception(__FILE__ , __LINE__ );
70  }
71 
72  PCAObj = std::make_unique<PCAHandler>();
73  //Check whether first and last pcadpar are set and if not just PCA everything
74  if(FirstPCAdpar == -999 || LastPCAdpar == -999){
75  if(FirstPCAdpar == -999 && LastPCAdpar == -999){
76  FirstPCAdpar = 0;
77  LastPCAdpar = covMatrix->GetNrows()-1;
78  }
79  else{
80  MACH3LOG_ERROR("You must either leave FirstPCAdpar and LastPCAdpar at -999 or set them both to something");
81  throw MaCh3Exception(__FILE__ , __LINE__ );
82  }
83  }
84 
85  PCAObj->ConstructPCA(covMatrix, FirstPCAdpar, LastPCAdpar, eigen_threshold, _fNumPar);
86  PCAObj->SetupPointers(&_fCurrVal, &_fPropVal);
87  // Make a note that we have now done PCA
88  pca = true;
89 }
90 
91 // ********************************************
92 void ParameterHandlerBase::Init(const std::string& name, const std::string& file) {
93 // ********************************************
94  // Set the covariance matrix from input ROOT file (e.g. flux, ND280, NIWG)
95  TFile *infile = new TFile(file.c_str(), "READ");
96  if (infile->IsZombie()) {
97  MACH3LOG_ERROR("Could not open input covariance ROOT file {} !!!", file);
98  MACH3LOG_ERROR("Was about to retrieve matrix with name {}", name);
99  throw MaCh3Exception(__FILE__ , __LINE__ );
100  }
101 
102  TMatrixDSym *CovMat = static_cast<TMatrixDSym*>(infile->Get(name.c_str()));
103 
104  if (!CovMat) {
105  MACH3LOG_ERROR("Could not find covariance matrix name {} in file {}", name, file);
106  MACH3LOG_ERROR("Are you really sure {} exists in the file?", name);
107  throw MaCh3Exception(__FILE__ , __LINE__ );
108  }
109 
110  PrintLength = 35;
111 
112  const int nThreads = M3::GetNThreads();
113  //KS: set Random numbers for each thread so each thread has different seed
114  //or for one thread if without MULTITHREAD
115  random_number.reserve(nThreads);
116  for (int iThread = 0; iThread < nThreads; iThread++) {
117  random_number.emplace_back(std::make_unique<TRandom3>(0));
118  }
119  // Not using adaptive by default
120  use_adaptive = false;
121  // Set the covariance matrix
122  _fNumPar = CovMat->GetNrows();
123 
124  InvertCovMatrix.resize(_fNumPar, std::vector<double>(_fNumPar, 0.0));
125  throwMatrixCholDecomp = new double*[_fNumPar]();
126  // Set the defaults to true
127  for(int i = 0; i < _fNumPar; i++) {
128  throwMatrixCholDecomp[i] = new double[_fNumPar]();
129  for (int j = 0; j < _fNumPar; j++) {
130  throwMatrixCholDecomp[i][j] = 0.;
131  }
132  }
133  SetName(name);
134  MakePosDef(CovMat);
135  SetCovMatrix(CovMat);
136  if (_fNumPar <= 0) {
137  MACH3LOG_CRITICAL("Covariance matrix {} has {} entries!", GetName(), _fNumPar);
138  throw MaCh3Exception(__FILE__ , __LINE__ );
139  }
140 
142 
143  infile->Close();
144 
145  MACH3LOG_INFO("Created covariance matrix named: {}", GetName());
146  MACH3LOG_INFO("from file: {}", file);
147  delete infile;
148 }
149 
150 // ********************************************
151 // ETA An init function for the YAML constructor
152 // All you really need from the YAML file is the number of Systematics
153 void ParameterHandlerBase::Init(const std::vector<std::string>& YAMLFile) {
154 // ********************************************
155  _fYAMLDoc["Systematics"] = YAML::Node(YAML::NodeType::Sequence);
156  for(unsigned int i = 0; i < YAMLFile.size(); i++)
157  {
158  YAML::Node YAMLDocTemp = M3OpenConfig(YAMLFile[i]);
159  for (const auto& item : YAMLDocTemp["Systematics"]) {
160  _fYAMLDoc["Systematics"].push_back(item);
161  }
162  }
163 
164  const int nThreads = M3::GetNThreads();
165  //KS: set Random numbers for each thread so each thread has different seed
166  //or for one thread if without MULTITHREAD
167  random_number.reserve(nThreads);
168  for (int iThread = 0; iThread < nThreads; iThread++) {
169  random_number.emplace_back(std::make_unique<TRandom3>(0));
170  }
171  PrintLength = 35;
172 
173  // Set the covariance matrix
174  _fNumPar = int(_fYAMLDoc["Systematics"].size());
175 
176  use_adaptive = false;
177 
178  InvertCovMatrix.resize(_fNumPar, std::vector<double>(_fNumPar, 0.0));
179  throwMatrixCholDecomp = new double*[_fNumPar]();
180  for(int i = 0; i < _fNumPar; i++) {
181  throwMatrixCholDecomp[i] = new double[_fNumPar]();
182  for (int j = 0; j < _fNumPar; j++) {
183  throwMatrixCholDecomp[i][j] = 0.;
184  }
185  }
187 
188  TMatrixDSym* _fCovMatrix = new TMatrixDSym(_fNumPar);
189  int i = 0;
190  std::vector<std::map<std::string,double>> Correlations(_fNumPar);
191  std::map<std::string, int> CorrNamesMap;
192 
193  //ETA - read in the systematics. Would be good to add in some checks to make sure
194  //that there are the correct number of entries i.e. are the _fNumPar for Names,
195  //PreFitValues etc etc.
196  for (auto const &param : _fYAMLDoc["Systematics"])
197  {
198  _fFancyNames[i] = Get<std::string>(param["Systematic"]["Names"]["FancyName"], __FILE__ , __LINE__);
199  _fPreFitValue[i] = Get<double>(param["Systematic"]["ParameterValues"]["PreFitValue"], __FILE__ , __LINE__);
200  _fIndivStepScale[i] = Get<double>(param["Systematic"]["StepScale"]["MCMC"], __FILE__ , __LINE__);
201  _fError[i] = Get<double>(param["Systematic"]["Error"], __FILE__ , __LINE__);
202  _fSampleNames[i] = GetFromManager<std::vector<std::string>>(param["Systematic"]["SampleNames"], {}, __FILE__, __LINE__);
203  if(_fError[i] <= 0) {
204  MACH3LOG_ERROR("Error for param {}({}) is negative and equal to {}", _fFancyNames[i], i, _fError[i]);
205  throw MaCh3Exception(__FILE__ , __LINE__ );
206  }
207  //ETA - a bit of a fudge but works
208  auto TempBoundsVec = GetBounds(param["Systematic"]["ParameterBounds"]);
209  _fLowBound[i] = TempBoundsVec[0];
210  _fUpBound[i] = TempBoundsVec[1];
211 
212  //ETA - now for parameters which are optional and have default values
213  _fFlatPrior[i] = GetFromManager<bool>(param["Systematic"]["FlatPrior"], false, __FILE__ , __LINE__);
214 
215  // 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
216  if(GetFromManager<bool>(param["Systematic"]["FixParam"], false, __FILE__ , __LINE__)) {
218  }
219 
220  if(param["Systematic"]["SpecialProposal"]) {
221  EnableSpecialProposal(param["Systematic"]["SpecialProposal"], i);
222  }
223 
224  //Fill the map to get the correlations later as well
225  CorrNamesMap[param["Systematic"]["Names"]["FancyName"].as<std::string>()]=i;
226 
227  //Also loop through the correlations
228  if(param["Systematic"]["Correlations"]) {
229  for(unsigned int Corr_i = 0; Corr_i < param["Systematic"]["Correlations"].size(); ++Corr_i){
230  for (YAML::const_iterator it = param["Systematic"]["Correlations"][Corr_i].begin(); it!=param["Systematic"]["Correlations"][Corr_i].end();++it) {
231  Correlations[i][it->first.as<std::string>()] = it->second.as<double>();
232  }
233  }
234  }
235  i++;
236  } // end loop over para
237  if(i != _fNumPar) {
238  MACH3LOG_CRITICAL("Inconsistent number of params in Yaml {} vs {}, this indicate wrong syntax", i, i, _fNumPar);
239  throw MaCh3Exception(__FILE__ , __LINE__ );
240  }
241  // ETA Now that we've been through all systematic let's fill the covmatrix
242  //This makes the root TCov from YAML
243  for(int j = 0; j < _fNumPar; j++) {
244  (*_fCovMatrix)(j, j) = _fError[j]*_fError[j];
245  //Get the map of parameter name to correlation from the Correlations object
246  for (auto const& pair : Correlations[j]) {
247  auto const& key = pair.first;
248  auto const& val = pair.second;
249  int index = -1;
250  //If you found the parameter name then get the index
251  if (CorrNamesMap.find(key) != CorrNamesMap.end()) {
252  index = CorrNamesMap[key];
253  } else {
254  MACH3LOG_ERROR("Parameter {} not in list! Check your spelling?", key);
255  throw MaCh3Exception(__FILE__ , __LINE__ );
256  }
257  double Corr1 = val;
258  double Corr2 = 0;
259  if(Correlations[index].find(_fFancyNames[j]) != Correlations[index].end()) {
260  Corr2 = Correlations[index][_fFancyNames[j]];
261  //Do they agree to better than float precision?
262  if(std::abs(Corr2 - Corr1) > FLT_EPSILON) {
263  MACH3LOG_ERROR("Correlations are not equal between {} and {}", _fFancyNames[j], key);
264  MACH3LOG_ERROR("Got : {} and {}", Corr2, Corr1);
265  throw MaCh3Exception(__FILE__ , __LINE__ );
266  }
267  } else {
268  MACH3LOG_ERROR("Correlation does not appear reciprocally between {} and {}", _fFancyNames[j], key);
269  throw MaCh3Exception(__FILE__ , __LINE__ );
270  }
271  (*_fCovMatrix)(j, index)= (*_fCovMatrix)(index, j) = Corr1*_fError[j]*_fError[index];
272  }
273  }
274 
275  //Now make positive definite
276  MakePosDef(_fCovMatrix);
277  SetCovMatrix(_fCovMatrix);
278 
279  if (_fNumPar <= 0) {
280  MACH3LOG_ERROR("ParameterHandler object has {} systematics!", _fNumPar);
281  throw MaCh3Exception(__FILE__ , __LINE__ );
282  }
283 
284  Tunes = std::make_unique<ParameterTunes>(_fYAMLDoc["Systematics"]);
285 
286  MACH3LOG_INFO("Created covariance matrix from files: ");
287  for(const auto &file : YAMLFile){
288  MACH3LOG_INFO("{} ", file);
289  }
290  MACH3LOG_INFO("----------------");
291  MACH3LOG_INFO("Found {} systematics parameters in total", _fNumPar);
292  MACH3LOG_INFO("----------------");
293 }
294 
295 // ********************************************
296 void ParameterHandlerBase::EnableSpecialProposal(const YAML::Node& param, const int Index){
297 // ********************************************
298  doSpecialStepProposal = true;
299 
300  bool CircEnabled = false;
301  bool FlipEnabled = false;
302 
303  if (param["CircularBounds"]) {
304  CircEnabled = true;
305  }
306 
307  if (param["FlipParameter"]) {
308  FlipEnabled = true;
309  }
310 
311  if (!CircEnabled && !FlipEnabled) {
312  MACH3LOG_ERROR("None of Special Proposal were enabled even though param {}, has SpecialProposal entry in Yaml", GetParFancyName(Index));
313  throw MaCh3Exception(__FILE__, __LINE__);
314  }
315 
316  if (CircEnabled) {
317  CircularBoundsIndex.push_back(Index);
318  CircularBoundsValues.push_back(Get<std::pair<double, double>>(param["CircularBounds"], __FILE__, __LINE__));
319  MACH3LOG_INFO("Enabling CircularBounds for parameter {} with range [{}, {}]",
320  GetParFancyName(Index),
321  CircularBoundsValues.back().first,
322  CircularBoundsValues.back().second);
323  // KS: Make sure circular bounds are within physical bounds. If we are outside of physics bound MCMC will never explore such phase space region
324  if (CircularBoundsValues.back().first < _fLowBound.at(Index) || CircularBoundsValues.back().second > _fUpBound.at(Index)) {
325  MACH3LOG_ERROR("Circular bounds [{}, {}] for parameter {} exceed physical bounds [{}, {}]",
326  CircularBoundsValues.back().first, CircularBoundsValues.back().second,
327  GetParFancyName(Index),
328  _fLowBound.at(Index), _fUpBound.at(Index));
329  throw MaCh3Exception(__FILE__, __LINE__);
330  }
331  }
332 
333  if (FlipEnabled) {
334  FlipParameterIndex.push_back(Index);
335  FlipParameterPoint.push_back(Get<double>(param["FlipParameter"], __FILE__, __LINE__));
336  MACH3LOG_INFO("Enabling Flipping for parameter {} with value {}",
337  GetParFancyName(Index),
338  FlipParameterPoint.back());
339  }
340 
341  if (CircEnabled && FlipEnabled) {
342  if (FlipParameterPoint.back() < CircularBoundsValues.back().first || FlipParameterPoint.back() > CircularBoundsValues.back().second) {
343  MACH3LOG_ERROR("FlipParameter value {} for parameter {} is outside the CircularBounds [{}, {}]",
344  FlipParameterPoint.back(), GetParFancyName(Index), CircularBoundsValues.back().first, CircularBoundsValues.back().second);
345  throw MaCh3Exception(__FILE__, __LINE__);
346  }
347 
348  const double low = CircularBoundsValues.back().first;
349  const double high = CircularBoundsValues.back().second;
350 
351  // Sanity check: ensure flipping any x in [low, high] keeps the result in [low, high]
352  const double flipped_low = 2 * FlipParameterPoint.back() - low;
353  const double flipped_high = 2 * FlipParameterPoint.back() - high;
354  const double min_flip = std::min(flipped_low, flipped_high);
355  const double max_flip = std::max(flipped_low, flipped_high);
356 
357  if (min_flip < low || max_flip > high) {
358  MACH3LOG_ERROR("Flipping about point {} for parameter {} would leave circular bounds [{}, {}]",
359  FlipParameterPoint.back(), GetParFancyName(Index), low, high);
360  throw MaCh3Exception(__FILE__, __LINE__);
361  }
362  }
363 }
364 
365 // ********************************************
366 // Set the covariance matrix for this class
367 void ParameterHandlerBase::SetCovMatrix(TMatrixDSym *cov) {
368 // ********************************************
369  if (cov == nullptr) {
370  MACH3LOG_ERROR("Could not find covariance matrix you provided to {}", __func__ );
371  throw MaCh3Exception(__FILE__ , __LINE__ );
372  }
373  covMatrix = cov;
374 
375  invCovMatrix = static_cast<TMatrixDSym *>(cov->Clone());
376  invCovMatrix->Invert();
377  //KS: ROOT has bad memory management, using standard double means we can decrease most operation by factor 2 simply due to cache hits
378  for (int i = 0; i < _fNumPar; i++)
379  {
380  for (int j = 0; j < _fNumPar; ++j)
381  {
382  InvertCovMatrix[i][j] = (*invCovMatrix)(i,j);
383  }
384  }
385 
386  SetThrowMatrix(cov);
387 }
388 // ********************************************
389 void ParameterHandlerBase::ReserveMemory(const int SizeVec) {
390 // ********************************************
391  _fNames = std::vector<std::string>(SizeVec);
392  _fFancyNames = std::vector<std::string>(SizeVec);
393  _fPreFitValue = std::vector<double>(SizeVec);
394  _fError = std::vector<double>(SizeVec);
395  _fCurrVal = std::vector<double>(SizeVec);
396  _fPropVal = std::vector<double>(SizeVec);
397  _fLowBound = std::vector<double>(SizeVec);
398  _fUpBound = std::vector<double>(SizeVec);
399  _fFlatPrior = std::vector<bool>(SizeVec);
400  _fIndivStepScale = std::vector<double>(SizeVec);
401  _fSampleNames = std::vector<std::vector<std::string>>(_fNumPar);
402 
403  corr_throw = new double[SizeVec];
404  // set random parameter vector (for correlated steps)
405  randParams = new double[SizeVec];
406 
407  // Set the defaults to true
408  for(int i = 0; i < SizeVec; i++) {
409  _fPreFitValue.at(i) = 1.;
410  _fError.at(i) = 1.;
411  _fCurrVal.at(i) = 0.;
412  _fPropVal.at(i) = 0.;
413  _fLowBound.at(i) = -999.99;
414  _fUpBound.at(i) = 999.99;
415  _fFlatPrior.at(i) = false;
416  _fIndivStepScale.at(i) = 1.;
417  corr_throw[i] = 0.0;
418  randParams[i] = 0.0;
419  }
420 
421  _fGlobalStepScale = 1.0;
422 }
423 
424 // ********************************************
425 // Set all the covariance matrix parameters to a user-defined value
426 // Might want to split this
427 void ParameterHandlerBase::SetPar(int i , double val) {
428 // ********************************************
429  MACH3LOG_INFO("Over-riding {}: ", GetParName(i));
430  MACH3LOG_INFO("_fPropVal ({}), _fCurrVal ({}), _fPreFitValue ({}) to ({})", _fPropVal[i], _fCurrVal[i], _fPreFitValue[i], val);
431 
432  _fPropVal[i] = val;
433  _fCurrVal[i] = val;
434  _fPreFitValue[i] = val;
435 
436  // Transfer the parameter values to the PCA basis
437  if (pca) PCAObj->TransferToPCA();
438 }
439 
440 // ********************************************
441 std::vector<double> ParameterHandlerBase::GetProposed() const {
442 // ********************************************
443  std::vector<double> props(_fNumPar);
444  for (int i = 0; i < _fNumPar; ++i) props[i] = _fPropVal[i];
445  return props;
446 }
447 
448 // *************************************
449 // Throw the parameters according to the covariance matrix
450 // This shouldn't be used in MCMC code ase it can break Detailed Balance;
452 // *************************************
453  // First draw new randParams
454  Randomize();
455 
457 
458  // KS: We use PCA very rarely on top PCA functionality isn't implemented for this function.
459  // Use __builtin_expect to give compiler a hint which option is more likely, which should help
460  // with better optimisation. This isn't critical but more to have example
461  if (__builtin_expect(!pca, 1)) {
462  #ifdef MULTITHREAD
463  #pragma omp parallel for
464  #endif
465  for (int i = 0; i < _fNumPar; ++i) {
466  // Check if parameter is fixed first: if so don't randomly throw
467  if (IsParameterFixed(i)) continue;
468 
469  _fPropVal[i] = _fPreFitValue[i] + corr_throw[i];
470  int throws = 0;
471  // Try again if we the initial parameter proposal falls outside of the range of the parameter
472  while (_fPropVal[i] > _fUpBound[i] || _fPropVal[i] < _fLowBound[i]) {
473  randParams[i] = random_number[M3::GetThreadIndex()]->Gaus(0, 1);
474  const double corr_throw_single = M3::MatrixVectorMultiSingle(throwMatrixCholDecomp, randParams, _fNumPar, i);
475  _fPropVal[i] = _fPreFitValue[i] + corr_throw_single;
476  if (throws > 10000)
477  {
478  //KS: Since we are multithreading there is danger that those messages
479  //will be all over the place, small price to pay for faster code
480  MACH3LOG_WARN("Tried {} times to throw parameter {} but failed", throws, i);
481  MACH3LOG_WARN("Matrix: {}", matrixName);
482  MACH3LOG_WARN("Param: {}", _fNames[i]);
483  MACH3LOG_WARN("Setting _fPropVal: {} to {}", _fPropVal[i], _fPreFitValue[i]);
484  MACH3LOG_WARN("I live at {}:{}", __FILE__, __LINE__);
485  _fPropVal[i] = _fPreFitValue[i];
486  //throw MaCh3Exception(__FILE__ , __LINE__ );
487  }
488  throws++;
489  }
490  _fCurrVal[i] = _fPropVal[i];
491  }
492  }
493  else
494  {
495  PCAObj->ThrowParameters(random_number, throwMatrixCholDecomp,
498  } // end if pca
499 }
500 
501 // *************************************
502 // Throw each parameter within their 1 sigma range
503 // Used to start the chain in different states
505 // *************************************
506  // Have the 1 sigma for each parameter in each covariance class, sweet!
507  // Don't want to change the prior array because that's what determines our likelihood
508  // Want to change the _fPropVal, _fCurrVal, _fPreFitValue
509  // _fPreFitValue and the others will already be set
510  for (int i = 0; i < _fNumPar; ++i) {
511  // Check if parameter is fixed first: if so don't randomly throw
512  if (IsParameterFixed(i)) continue;
513  // Check that the sigma range is larger than the parameter range
514  // If not, throw in the valid parameter range instead
515  const double paramrange = _fUpBound[i] - _fLowBound[i];
516  const double sigma = sqrt((*covMatrix)(i,i));
517  double throwrange = sigma;
518  if (paramrange < sigma) throwrange = paramrange;
519 
520  _fPropVal[i] = _fPreFitValue[i] + random_number[0]->Gaus(0, 1)*throwrange;
521  // Try again if we the initial parameter proposal falls outside of the range of the parameter
522  int throws = 0;
523  while (_fPropVal[i] > _fUpBound[i] || _fPropVal[i] < _fLowBound[i]) {
524  if (throws > 1000) {
525  MACH3LOG_WARN("Tried {} times to throw parameter {} but failed", throws, i);
526  MACH3LOG_WARN("Matrix: {}", matrixName);
527  MACH3LOG_WARN("Param: {}", _fNames[i]);
528  throw MaCh3Exception(__FILE__ , __LINE__ );
529  }
530  _fPropVal[i] = _fPreFitValue[i] + random_number[0]->Gaus(0, 1)*throwrange;
531  throws++;
532  }
533  MACH3LOG_INFO("Setting current step in {} param {} = {} from {}", matrixName, i, _fPropVal[i], _fCurrVal[i]);
534  _fCurrVal[i] = _fPropVal[i];
535  }
536  if (pca) PCAObj->TransferToPCA();
537 }
538 
539 // *************************************
540 // Set a single parameter
541 void ParameterHandlerBase::SetSingleParameter(const int parNo, const double parVal) {
542 // *************************************
543  _fPropVal[parNo] = parVal;
544  _fCurrVal[parNo] = parVal;
545  MACH3LOG_DEBUG("Setting {} (parameter {}) to {})", GetParName(parNo), parNo, parVal);
546  if (pca) PCAObj->TransferToPCA();
547 }
548 
549 // ********************************************
550 void ParameterHandlerBase::SetParCurrProp(const int parNo, const double parVal) {
551 // ********************************************
552  _fPropVal[parNo] = parVal;
553  _fCurrVal[parNo] = parVal;
554  MACH3LOG_DEBUG("Setting {} (parameter {}) to {})", GetParName(parNo), parNo, parVal);
555  if (pca) PCAObj->TransferToPCA();
556 }
557 
558 // ************************************************
559 // Propose a step for the set of systematics parameters this covariance class holds
561 // ************************************************
562  // Make the random numbers for the step proposal
563  Randomize();
564  CorrelateSteps();
565 
566  // KS: According to Dr Wallace we update using previous not proposed step
567  // this way we do special proposal after adaptive after.
568  // This way we can shortcut and skip rest of proposal
569  if(!doSpecialStepProposal) return;
570 
572 }
573 
574 // ************************************************
576 // ************************************************
578 
579  // HW It should now automatically set dcp to be with [-pi, pi]
580  for (size_t i = 0; i < CircularBoundsIndex.size(); ++i) {
581  const int index = CircularBoundsIndex[i];
582  if(!IsParameterFixed(index))
583  CircularParBounds(index, CircularBoundsValues[i].first, CircularBoundsValues[i].second);
584  }
585 
586  // Okay now we've done the standard steps, we can add in our nice flips hierarchy flip first
587  for (size_t i = 0; i < FlipParameterIndex.size(); ++i) {
588  const int index = FlipParameterIndex[i];
589  if(!IsParameterFixed(index))
591  }
592 }
593 
594 // ************************************************
595 // "Randomize" the parameters in the covariance class for the proposed step
596 // Used the proposal kernel and the current parameter value to set proposed step
597 // Also get a new random number for the randParams
599 // ************************************************
600  if (!pca) {
601  //KS: By multithreading here we gain at least factor 2 with 8 threads with ND only fit
602  #ifdef MULTITHREAD
603  #pragma omp parallel for
604  #endif
605  for (int i = 0; i < _fNumPar; ++i) {
606  // If parameter isn't fixed
607  if (!IsParameterFixed(i) > 0.0) {
608  randParams[i] = random_number[M3::GetThreadIndex()]->Gaus(0, 1);
609  // If parameter IS fixed
610  } else {
611  randParams[i] = 0.0;
612  }
613  } // end for
614  // If we're in the PCA basis we instead throw parameters there (only _fNumParPCA parameter)
615  } else {
616  // Scale the random parameters by the sqrt of eigen values for the throw
617  #ifdef MULTITHREAD
618  #pragma omp parallel for
619  #endif
620  for (int i = 0; i < PCAObj->GetNumberPCAedParameters(); ++i)
621  {
622  // If parameter IS fixed or out of bounds
623  if (PCAObj->IsParameterFixedPCA(i)) {
624  randParams[i] = 0.0;
625  } else {
626  randParams[i] = random_number[M3::GetThreadIndex()]->Gaus(0,1);
627  }
628  }
629  }
630 }
631 
632 // ************************************************
633 // Correlate the steps by setting the proposed step of a parameter to its current value + some correlated throw
635 // ************************************************
636  //KS: Using custom function compared to ROOT one with 8 threads we have almost factor 2 performance increase, by replacing TMatrix with just double we increase it even more
638 
639  // If not doing PCA
640  if (!pca) {
641  #ifdef MULTITHREAD
642  #pragma omp parallel for
643  #endif
644  for (int i = 0; i < _fNumPar; ++i) {
645  if (!IsParameterFixed(i) > 0.) {
647  }
648  }
649  // If doing PCA throw uncorrelated in PCA basis (orthogonal basis by definition)
650  } else {
652  }
653 }
654 // ********************************************
655 // Update so that current step becomes the previously proposed step
657 // ********************************************
658  if (!pca) {
659  #ifdef MULTITHREAD
660  #pragma omp parallel for
661  #endif
662  for (int i = 0; i < _fNumPar; ++i) {
663  // Update state so that current state is proposed state
664  _fCurrVal[i] = _fPropVal[i];
665  }
666  } else {
667  PCAObj->AcceptStep();
668  }
669 }
670 
671 // *************************************
672 //HW: This method is a tad hacky but modular arithmetic gives me a headache.
673 void ParameterHandlerBase::CircularParBounds(const int index, const double LowBound, const double UpBound) {
674 // *************************************
675  if(_fPropVal[index] > UpBound) {
676  _fPropVal[index] = LowBound + std::fmod(_fPropVal[index] - UpBound, UpBound - LowBound);
677  } else if (_fPropVal[index] < LowBound) {
678  _fPropVal[index] = UpBound - std::fmod(LowBound - _fPropVal[index], UpBound - LowBound);
679  }
680 }
681 
682 // *************************************
683 void ParameterHandlerBase::FlipParameterValue(const int index, const double FlipPoint) {
684 // *************************************
685  if(random_number[0]->Uniform() < 0.5) {
686  _fPropVal[index] = 2 * FlipPoint - _fPropVal[index];
687  }
688 }
689 
690 // ********************************************
691 // Throw the proposed parameter by mag sigma
692 // Should really just have the user specify this throw by having argument double
693 void ParameterHandlerBase::ThrowParProp(const double mag) {
694 // ********************************************
695  Randomize();
696  if (!pca) {
697  // Make the correlated throw
699  // Number of sigmas we throw
700  for (int i = 0; i < _fNumPar; i++) {
701  if (!IsParameterFixed(i) > 0.)
702  _fPropVal[i] = _fCurrVal[i] + corr_throw[i]*mag;
703  }
704  } else {
705  PCAObj->ThrowParProp(mag, randParams);
706  }
707 }
708 // ********************************************
709 // Helper function to throw the current parameter by mag sigmas
710 // Can study bias in MCMC with this; put different starting parameters
711 void ParameterHandlerBase::ThrowParCurr(const double mag) {
712 // ********************************************
713  Randomize();
714  if (!pca) {
715  // Get the correlated throw vector
717  // The number of sigmas to throw
718  // Should probably have this as a default parameter input to the function instead
719  for (int i = 0; i < _fNumPar; i++) {
720  if (!IsParameterFixed(i) > 0.){
721  _fCurrVal[i] = corr_throw[i]*mag;
722  }
723  }
724  } else {
725  PCAObj->ThrowParCurr(mag, randParams);
726  }
727 }
728 // ********************************************
729 // Function to print the prior values
731 // ********************************************
732  MACH3LOG_INFO("Prior values for {} ParameterHandler:", GetName());
733  for (int i = 0; i < _fNumPar; i++) {
734  MACH3LOG_INFO(" {} {} ", GetParFancyName(i), GetParInit(i));
735  }
736 }
737 
738 // ********************************************
739 // Function to print the prior, current and proposed values
741 // ********************************************
742  MACH3LOG_INFO("Printing parameters for {}", GetName());
743  // Dump out the PCA parameters too
744  if (pca) {
745  PCAObj->Print();
746  }
747  MACH3LOG_INFO("{:<30} {:<10} {:<10} {:<10}", "Name", "Prior", "Current", "Proposed");
748  for (int i = 0; i < _fNumPar; ++i) {
749  MACH3LOG_INFO("{:<30} {:<10.2f} {:<10.2f} {:<10.2f}", GetParFancyName(i), _fPreFitValue[i], _fCurrVal[i], _fPropVal[i]);
750  }
751 }
752 
753 // ********************************************
754 // Get the likelihood in the case where we want to include priors on the parameters
755 // fParEvalLikelihood stores if we want to evaluate the likelihood for the given parameter
756 // true = evaluate likelihood (so run with a prior)
757 // false = don't evaluate likelihood (so run without a prior)
759 // ********************************************
760  double logL = 0.0;
761  #ifdef MULTITHREAD
762  #pragma omp parallel for reduction(+:logL)
763  #endif
764  for(int i = 0; i < _fNumPar; ++i){
765  if(_fFlatPrior[i]){
766  //HW: Flat prior, no need to calculate anything
767  continue;
768  }
769 
770  #ifdef MULTITHREAD
771  #pragma omp simd
772  #endif
773  for (int j = 0; j <= i; ++j) {
774  if (!_fFlatPrior[j]) {
775  //KS: Since matrix is symmetric we can calculate non diagonal elements only once and multiply by 2, can bring up to factor speed decrease.
776  double scale = (i != j) ? 1. : 0.5;
777  logL += scale * (_fPropVal[i] - _fPreFitValue[i])*(_fPropVal[j] - _fPreFitValue[j])*InvertCovMatrix[i][j];
778  }
779  }
780  }
781  return logL;
782 }
783 
784 // ********************************************
786 // ********************************************
787  int NOutside = 0;
788  #ifdef MULTITHREAD
789  #pragma omp parallel for reduction(+:NOutside)
790  #endif
791  for (int i = 0; i < _fNumPar; ++i){
792  if(_fPropVal[i] > _fUpBound[i] || _fPropVal[i] < _fLowBound[i]){
793  NOutside++;
794  }
795  }
796  return NOutside;
797 }
798 
799 // ********************************************
801 // ********************************************
802  // Default behaviour is to reject negative values + do std llh calculation
803  const int NOutside = CheckBounds();
804 
805  if(NOutside > 0) return NOutside*M3::_LARGE_LOGL_;
806 
807  return CalcLikelihood();
808 }
809 
810 // ********************************************
811 // Sets the proposed parameters to the prior values
812 void ParameterHandlerBase::SetParameters(const std::vector<double>& pars) {
813 // ********************************************
814  // If empty, set the proposed to prior
815  if (pars.empty()) {
816  // For xsec this means setting to the prior (because prior is the prior)
817  for (int i = 0; i < _fNumPar; i++) {
818  _fPropVal[i] = _fPreFitValue[i];
819  }
820  // If not empty, set the parameters to the specified
821  } else {
822  if (pars.size() != size_t(_fNumPar)) {
823  MACH3LOG_ERROR("Parameter arrays of incompatible size! Not changing parameters! {} has size {} but was expecting {}", matrixName, pars.size(), _fNumPar);
824  throw MaCh3Exception(__FILE__ , __LINE__ );
825  }
826  int parsSize = int(pars.size());
827  for (int i = 0; i < parsSize; i++) {
828  //Make sure that you are actually passing a number to set the parameter to
829  if(std::isnan(pars[i])) {
830  MACH3LOG_ERROR("Trying to set parameter value to a nan for parameter {} in matrix {}. This will not go well!", GetParName(i), matrixName);
831  throw MaCh3Exception(__FILE__ , __LINE__ );
832  } else {
833  _fPropVal[i] = pars[i];
834  }
835  }
836  }
837  // And if pca make the transfer
838  if (pca) {
839  PCAObj->TransferToPCA();
840  PCAObj->TransferToParam();
841  }
842 }
843 
844 // ********************************************
845 void ParameterHandlerBase::SetBranches(TTree &tree, bool SaveProposal) {
846 // ********************************************
847  // loop over parameters and set a branch
848  for (int i = 0; i < _fNumPar; ++i) {
849  tree.Branch(_fNames[i].c_str(), &_fCurrVal[i], Form("%s/D", _fNames[i].c_str()));
850  }
851  // When running PCA, also save PCA parameters
852  if (pca) {
853  PCAObj->SetBranches(tree, SaveProposal, _fNames);
854  }
855  if(SaveProposal)
856  {
857  // loop over parameters and set a branch
858  for (int i = 0; i < _fNumPar; ++i) {
859  tree.Branch(Form("%s_Prop", _fNames[i].c_str()), &_fPropVal[i], Form("%s_Prop/D", _fNames[i].c_str()));
860  }
861  }
862 }
863 
864 // ********************************************
865 void ParameterHandlerBase::SetStepScale(const double scale, const bool verbose) {
866 // ********************************************
867  if(scale <= 0) {
868  MACH3LOG_ERROR("You are trying so set StepScale to 0 or negative this will not work");
869  throw MaCh3Exception(__FILE__ , __LINE__ );
870  }
871 
872  if(verbose){
873  MACH3LOG_INFO("{} setStepScale() = {}", GetName(), scale);
874  const double SuggestedScale = 2.38*2.38/_fNumPar;
875  if(std::fabs(scale - SuggestedScale)/SuggestedScale > 1) {
876  MACH3LOG_WARN("Defined Global StepScale is {}, while suggested suggested {}", scale, SuggestedScale);
877  }
878  }
879  _fGlobalStepScale = scale;
880 }
881 
882 // ********************************************
883 int ParameterHandlerBase::GetParIndex(const std::string& name) const {
884 // ********************************************
885  int Index = M3::_BAD_INT_;
886  for (int i = 0; i <_fNumPar; ++i) {
887  if(name == _fFancyNames[i]) {
888  Index = i;
889  break;
890  }
891  }
892  return Index;
893 }
894 
895 // ********************************************
897 // ********************************************
898  // fix or unfix all parameters by multiplying by -1
899  if(!pca) {
900  for (int i = 0; i < _fNumPar; i++) _fError[i] *= -1.0;
901  } else{
902  PCAObj->ToggleFixAllParameters();
903  }
904 }
905 
906 // ********************************************
908 // ********************************************
909  if(!pca) {
910  if (i > _fNumPar) {
911  MACH3LOG_ERROR("Can't {} for parameter {} because size of covariance ={}", __func__, i, _fNumPar);
912  MACH3LOG_ERROR("Fix this in your config file please!");
913  throw MaCh3Exception(__FILE__ , __LINE__ );
914  } else {
915  _fError[i] *= -1.0;
916  MACH3LOG_INFO("Setting {}(parameter {}) to fixed at {}", GetParFancyName(i), i, _fCurrVal[i]);
917  }
918  } else {
919  PCAObj->ToggleFixParameter(i, _fNames);
920  }
921 }
922 
923 // ********************************************
924 void ParameterHandlerBase::ToggleFixParameter(const std::string& name) {
925 // ********************************************
926  const int Index = GetParIndex(name);
927  if(Index != M3::_BAD_INT_) {
928  ToggleFixParameter(Index);
929  return;
930  }
931 
932  MACH3LOG_WARN("I couldn't find parameter with name {}, therefore will not fix it", name);
933 }
934 
935 // ********************************************
936 bool ParameterHandlerBase::IsParameterFixed(const std::string& name) const {
937 // ********************************************
938  const int Index = GetParIndex(name);
939  if(Index != M3::_BAD_INT_) {
940  return IsParameterFixed(Index);
941  }
942 
943  MACH3LOG_WARN("I couldn't find parameter with name {}, therefore don't know if it fixed", name);
944  return false;
945 }
946 
947 // ********************************************
948 void ParameterHandlerBase::SetFlatPrior(const int i, const bool eL) {
949 // ********************************************
950  if (i > _fNumPar) {
951  MACH3LOG_INFO("Can't {} for Cov={}/Param={} because size of Covariance = {}", __func__, GetName(), i, _fNumPar);
952  MACH3LOG_ERROR("Fix this in your config file please!");
953  throw MaCh3Exception(__FILE__ , __LINE__ );
954  } else {
955  if(eL){
956  MACH3LOG_INFO("Setting {} (parameter {}) to flat prior", GetParName(i), i);
957  }
958  else{
959  // HW :: This is useful
960  MACH3LOG_INFO("Setting {} (parameter {}) to non-flat prior", GetParName(i), i);
961  }
962  _fFlatPrior[i] = eL;
963  }
964 }
965 
966 // ********************************************
967 void ParameterHandlerBase::SetIndivStepScale(const std::vector<double>& stepscale) {
968 // ********************************************
969  if (int(stepscale.size()) != _fNumPar)
970  {
971  MACH3LOG_WARN("Stepscale vector not equal to number of parameters. Quitting..");
972  MACH3LOG_WARN("Size of argument vector: {}", stepscale.size());
973  MACH3LOG_WARN("Expected size: {}", _fNumPar);
974  return;
975  }
976 
977  for (int iParam = 0 ; iParam < _fNumPar; iParam++) {
978  _fIndivStepScale[iParam] = stepscale[iParam];
979  }
981 }
982 
983 // ********************************************
985 // ********************************************
986  MACH3LOG_INFO("============================================================");
987  MACH3LOG_INFO("{:<{}} | {:<11}", "Parameter:", PrintLength, "Step scale:");
988  for (int iParam = 0; iParam < _fNumPar; iParam++) {
989  MACH3LOG_INFO("{:<{}} | {:<11}", _fFancyNames[iParam].c_str(), PrintLength, _fIndivStepScale[iParam]);
990  }
991  MACH3LOG_INFO("============================================================");
992 }
993 
994 // ********************************************
995 //Makes sure that matrix is positive-definite by adding a small number to on-diagonal elements
996 void ParameterHandlerBase::MakePosDef(TMatrixDSym *cov) {
997 // ********************************************
998  if(cov == nullptr){
999  cov = &*covMatrix;
1000  MACH3LOG_WARN("Passed nullptr to cov matrix in {}", matrixName);
1001  }
1002 
1003  M3::MakeMatrixPosDef(cov);
1004 }
1005 
1006 // ********************************************
1008 // ********************************************
1009  std::vector<double> stepScales(_fNumPar);
1010  for (int i = 0; i <_fNumPar; i++) {
1011  stepScales[i] = 1.;
1012  }
1013  _fGlobalStepScale = 1.0;
1014  SetIndivStepScale(stepScales);
1015 }
1016 
1017 // ********************************************
1018 // HW: Code for throwing from separate throw matrix, needs to be set after init to ensure pos-def
1020 // ********************************************
1021  if (cov == nullptr) {
1022  MACH3LOG_ERROR("Could not find covariance matrix you provided to {}", __func__);
1023  throw MaCh3Exception(__FILE__ , __LINE__ );
1024  }
1025 
1026  if (covMatrix->GetNrows() != cov->GetNrows()) {
1027  MACH3LOG_ERROR("Matrix given for throw Matrix is not the same size as the covariance matrix stored in object!");
1028  MACH3LOG_ERROR("Stored covariance matrix size: {}", covMatrix->GetNrows());
1029  MACH3LOG_ERROR("Given matrix size: {}", cov->GetNrows());
1030  throw MaCh3Exception(__FILE__ , __LINE__ );
1031  }
1032 
1033  throwMatrix = static_cast<TMatrixDSym*>(cov->Clone());
1034  if(use_adaptive && AdaptiveHandler->AdaptionUpdate()) MakeClosestPosDef(throwMatrix);
1035  else MakePosDef(throwMatrix);
1036 
1037  auto throwMatrix_CholDecomp = M3::GetCholeskyDecomposedMatrix(*throwMatrix, matrixName);
1038 
1039  //KS: ROOT has bad memory management, using standard double means we can decrease most operation by factor 2 simply due to cache hits
1040  #ifdef MULTITHREAD
1041  #pragma omp parallel for collapse(2)
1042  #endif
1043  for (int i = 0; i < _fNumPar; ++i)
1044  {
1045  for (int j = 0; j < _fNumPar; ++j)
1046  {
1047  throwMatrixCholDecomp[i][j] = throwMatrix_CholDecomp[i][j];
1048  }
1049  }
1050 }
1051 
1052 // ********************************************
1054 // ********************************************
1055  delete throwMatrix;
1056  throwMatrix = nullptr;
1057  SetThrowMatrix(cov);
1058 }
1059 
1060 // ********************************************
1061 // HW : Here be adaption
1062 void ParameterHandlerBase::InitialiseAdaption(const YAML::Node& adapt_manager){
1063 // ********************************************
1064  if(PCAObj){
1065  MACH3LOG_ERROR("PCA has been enabled and now trying to enable Adaption. Right now both configuration don't work with each other");
1066  throw MaCh3Exception(__FILE__ , __LINE__ );
1067  }
1068  if(AdaptiveHandler){
1069  MACH3LOG_ERROR("Adaptive Handler has already been initialise can't do it again so skipping.");
1070  return;
1071  }
1072  AdaptiveHandler = std::make_unique<adaptive_mcmc::AdaptiveMCMCHandler>();
1073  // Now we read the general settings [these SHOULD be common across all matrices!]
1074  bool success = AdaptiveHandler->InitFromConfig(adapt_manager, matrixName, &_fCurrVal, &_fError);
1075  if(!success) return;
1076  AdaptiveHandler->Print();
1077 
1078  // Next let"s check for external matrices
1079  // We"re going to grab this info from the YAML manager
1080  if(!GetFromManager<bool>(adapt_manager["AdaptionOptions"]["Covariance"][matrixName]["UseExternalMatrix"], false, __FILE__ , __LINE__)) {
1081  MACH3LOG_WARN("Not using external matrix for {}, initialising adaption from scratch", matrixName);
1082  // If we don't have a covariance matrix to start from for adaptive tune we need to make one!
1083  use_adaptive = true;
1084  AdaptiveHandler->CreateNewAdaptiveCovariance();
1085  return;
1086  }
1087 
1088  // Finally, we accept that we want to read the matrix from a file!
1089  auto external_file_name = GetFromManager<std::string>(adapt_manager["AdaptionOptions"]["Covariance"][matrixName]["ExternalMatrixFileName"], "", __FILE__ , __LINE__);
1090  auto external_matrix_name = GetFromManager<std::string>(adapt_manager["AdaptionOptions"]["Covariance"][matrixName]["ExternalMatrixName"], "", __FILE__ , __LINE__);
1091  auto external_mean_name = GetFromManager<std::string>(adapt_manager["AdaptionOptions"]["Covariance"][matrixName]["ExternalMeansName"], "", __FILE__ , __LINE__);
1092 
1093  AdaptiveHandler->SetThrowMatrixFromFile(external_file_name, external_matrix_name, external_mean_name, use_adaptive);
1094  SetThrowMatrix(AdaptiveHandler->GetAdaptiveCovariance());
1096 
1097  MACH3LOG_INFO("Successfully Set External Throw Matrix Stored in {}", external_file_name);
1098 }
1099 
1100 // ********************************************
1101 // Truely adaptive MCMC!
1103 // ********************************************
1104  // Updates adaptive matrix
1105  // First we update the total means
1106 
1107  // Skip this if we're at a large number of steps
1108  if(AdaptiveHandler->SkipAdaption()) {
1109  AdaptiveHandler->IncrementNSteps();
1110  return;
1111  }
1112 
1113  // Call main adaption function
1114  AdaptiveHandler->UpdateAdaptiveCovariance();
1115 
1116  //This is likely going to be the slow bit!
1117  if(AdaptiveHandler->IndivStepScaleAdapt()) {
1119  }
1120 
1121  if(AdaptiveHandler->UpdateMatrixAdapt()) {
1122  TMatrixDSym* update_matrix = static_cast<TMatrixDSym*>(AdaptiveHandler->GetAdaptiveCovariance()->Clone());
1123  UpdateThrowMatrix(update_matrix); //Now we update and continue!
1124  //Also Save the adaptive to file
1125  AdaptiveHandler->SaveAdaptiveToFile(AdaptiveHandler->GetOutFileName(), GetName());
1126  }
1127 
1128  AdaptiveHandler->IncrementNSteps();
1129 }
1130 
1131 // ********************************************
1132 //HW: Finds closest possible positive definite matrix in Frobenius Norm ||.||_frob
1133 // Where ||X||_frob=sqrt[sum_ij(x_ij^2)] (basically just turns an n,n matrix into vector in n^2 space
1134 // then does Euclidean norm)
1136 // ********************************************
1137  // Want to get cov' = (cov_sym+cov_polar)/2
1138  // cov_sym=(cov+cov^T)/2
1139  // cov_polar-> SVD cov to cov=USV^T then cov_polar=VSV^T
1140 
1141  //Get frob norm of cov
1142  // Double_t cov_norm=cov->E2Norm();
1143 
1144  TMatrixDSym* cov_trans = cov;
1145  cov_trans->T();
1146  TMatrixDSym cov_sym = 0.5*(*cov+*cov_trans); //If cov is symmetric does nothing, otherwise just ensures symmetry
1147 
1148  //Do SVD to get polar form
1149  TDecompSVD cov_sym_svd=TDecompSVD(cov_sym);
1150  if(!cov_sym_svd.Decompose()){
1151  MACH3LOG_ERROR("Cannot do SVD on input matrix!");
1152  throw MaCh3Exception(__FILE__ , __LINE__ );
1153  }
1154 
1155  TMatrixD cov_sym_v = cov_sym_svd.GetV();
1156  TMatrixD cov_sym_vt = cov_sym_v;
1157  cov_sym_vt.T();
1158  //SVD returns as vector (grrr) so need to get into matrix form for multiplying!
1159  TVectorD cov_sym_sigvect = cov_sym_svd.GetSig();
1160 
1161  const Int_t nCols = cov_sym_v.GetNcols(); //square so only need rows hence lack of cols
1162  TMatrixDSym cov_sym_sig(nCols);
1163  TMatrixDDiag cov_sym_sig_diag(cov_sym_sig);
1164  cov_sym_sig_diag=cov_sym_sigvect;
1165 
1166  //Can finally get H=VSV
1167  TMatrixDSym cov_sym_polar = cov_sym_sig.SimilarityT(cov_sym_vt);//V*S*V^T (this took forver to find!)
1168 
1169  //Now we can construct closest approximater Ahat=0.5*(B+H)
1170  TMatrixDSym cov_closest_approx = 0.5*(cov_sym+cov_sym_polar);//Not fully sure why this is even needed since symmetric B -> U=V
1171  //Get norm of transformed
1172  // Double_t approx_norm=cov_closest_approx.E2Norm();
1173  //MACH3LOG_INFO("Initial Norm: {:.6f} | Norm after transformation: {:.6f} | Ratio: {:.6f}", cov_norm, approx_norm, cov_norm / approx_norm);
1174 
1175  *cov = cov_closest_approx;
1176  //Now can just add a makeposdef!
1177  MakePosDef(cov);
1178 }
1179 
1180 // ********************************************
1181 // KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plotting
1183 // ********************************************
1184  TH2D* hMatrix = new TH2D(GetName().c_str(), GetName().c_str(), _fNumPar, 0.0, _fNumPar, _fNumPar, 0.0, _fNumPar);
1185  hMatrix->SetDirectory(nullptr);
1186  for(int i = 0; i < _fNumPar; i++)
1187  {
1188  hMatrix->SetBinContent(i+1, i+1, 1.);
1189  hMatrix->GetXaxis()->SetBinLabel(i+1, GetParFancyName(i).c_str());
1190  hMatrix->GetYaxis()->SetBinLabel(i+1, GetParFancyName(i).c_str());
1191  }
1192 
1193  #ifdef MULTITHREAD
1194  #pragma omp parallel for
1195  #endif
1196  for(int i = 0; i < _fNumPar; i++)
1197  {
1198  for(int j = 0; j <= i; j++)
1199  {
1200  const double Corr = (*covMatrix)(i,j) / ( GetDiagonalError(i) * GetDiagonalError(j));
1201  hMatrix->SetBinContent(i+1, j+1, Corr);
1202  hMatrix->SetBinContent(j+1, i+1, Corr);
1203  }
1204  }
1205  return hMatrix;
1206 }
1207 
1208 // ********************************************
1209 // KS: After step scale, prefit etc. value were modified save this modified config.
1211 // ********************************************
1212  if (!_fYAMLDoc)
1213  {
1214  MACH3LOG_CRITICAL("Yaml node hasn't been initialised for matrix {}, something is not right", matrixName);
1215  MACH3LOG_CRITICAL("I am not throwing error but should be investigated");
1216  return;
1217  }
1218 
1219  YAML::Node copyNode = _fYAMLDoc;
1220  int i = 0;
1221 
1222  for (YAML::Node param : copyNode["Systematics"])
1223  {
1224  //KS: Feel free to update it, if you need updated prefit value etc
1225  param["Systematic"]["StepScale"]["MCMC"] = MaCh3Utils::FormatDouble(_fIndivStepScale[i], 4);
1226  i++;
1227  }
1228  // Save the modified node to a file
1229  std::ofstream fout("Modified_Matrix.yaml");
1230  fout << copyNode;
1231  fout.close();
1232 }
1233 
1234 // ********************************************
1235 bool ParameterHandlerBase::AppliesToSample(const int SystIndex, const std::string& SampleName) const {
1236 // ********************************************
1237  // Empty means apply to all
1238  if (_fSampleNames[SystIndex].size() == 0) return true;
1239 
1240  // Make a copy and to lower case to not be case sensitive
1241  std::string SampleNameCopy = SampleName;
1242  std::transform(SampleNameCopy.begin(), SampleNameCopy.end(), SampleNameCopy.begin(), ::tolower);
1243 
1244  // Check for unsupported wildcards in SampleNameCopy
1245  if (SampleNameCopy.find('*') != std::string::npos) {
1246  MACH3LOG_ERROR("Wildcards ('*') are not supported in sample name: '{}'", SampleName);
1247  throw MaCh3Exception(__FILE__ , __LINE__ );
1248  }
1249 
1250  bool Applies = false;
1251 
1252  for (size_t i = 0; i < _fSampleNames[SystIndex].size(); i++) {
1253  // Convert to low case to not be case sensitive
1254  std::string pattern = _fSampleNames[SystIndex][i];
1255  std::transform(pattern.begin(), pattern.end(), pattern.begin(), ::tolower);
1256 
1257  // Replace '*' in the pattern with '.*' for regex matching
1258  std::string regexPattern = "^" + std::regex_replace(pattern, std::regex("\\*"), ".*") + "$";
1259  try {
1260  std::regex regex(regexPattern);
1261  if (std::regex_match(SampleNameCopy, regex)) {
1262  Applies = true;
1263  break;
1264  }
1265  } catch (const std::regex_error& e) {
1266  // Handle regex error (for invalid patterns)
1267  MACH3LOG_ERROR("Regex error: {}", e.what());
1268  }
1269  }
1270  return Applies;
1271 }
1272 
1273 // ********************************************
1274 // Set proposed parameter values vector to be base on tune values
1275 void ParameterHandlerBase::SetTune(const std::string& TuneName) {
1276 // ********************************************
1277  if(Tunes == nullptr) {
1278  MACH3LOG_ERROR("Tunes haven't been initialised, which are being loaded from YAML, have you used some deprecated constructor");
1279  throw MaCh3Exception(__FILE__, __LINE__);
1280  }
1281  auto Values = Tunes->GetTune(TuneName);
1282 
1283  SetParameters(Values);
1284 }
1285 
1286 // *************************************
1296  std::vector<double>& BranchValues,
1297  std::vector<std::string>& BranchNames) {
1298 // *************************************
1299  BranchValues.resize(GetNumParams());
1300  BranchNames.resize(GetNumParams());
1301 
1302  for (int i = 0; i < GetNumParams(); ++i) {
1303  BranchNames[i] = GetParName(i);
1304  if (!PosteriorFile->GetBranch(BranchNames[i].c_str())) {
1305  MACH3LOG_ERROR("Branch '{}' does not exist in the TTree!", BranchNames[i]);
1306  throw MaCh3Exception(__FILE__, __LINE__);
1307  }
1308  PosteriorFile->SetBranchStatus(BranchNames[i].c_str(), true);
1309  PosteriorFile->SetBranchAddress(BranchNames[i].c_str(), &BranchValues[i]);
1310  }
1311 }
#define _noexcept_
KS: noexcept can help with performance but is terrible for debugging, this is meant to help easy way ...
Definition: Core.h:86
int size
std::vector< TString > BranchNames
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:28
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:24
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
Type Get(const YAML::Node &node, const std::string File, const int Line)
Get content of config file.
Definition: YamlHelper.h:266
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Definition: YamlHelper.h:561
#define GetBounds(filename)
Definition: YamlHelper.h:562
Custom exception class for MaCh3 errors.
std::vector< int > CircularBoundsIndex
Indices of parameters with circular bounds.
void ToggleFixParameter(const int i)
fix parameter at prior values
bool use_adaptive
Are we using AMCMC?
virtual void ProposeStep()
Generate a new proposed state.
double * randParams
Random number taken from gaussian around prior error used for corr_throw.
double ** throwMatrixCholDecomp
Throw matrix that is being used in the fit, much faster as TMatrixDSym cache miss.
std::unique_ptr< adaptive_mcmc::AdaptiveMCMCHandler > AdaptiveHandler
Struct containing information about adaption.
virtual ~ParameterHandlerBase()
Destructor.
TMatrixDSym * invCovMatrix
The inverse covariance matrix.
std::vector< int > FlipParameterIndex
Indices of parameters with flip symmetry.
bool AppliesToSample(const int SystIndex, const std::string &SampleName) const
Check if parameter is affecting given sample name.
void AcceptStep() _noexcept_
Accepted this step.
void ThrowParProp(const double mag=1.)
Throw the proposed parameter by mag sigma. Should really just have the user specify this throw by hav...
std::unique_ptr< ParameterTunes > Tunes
Struct containing information about adaption.
void InitialiseAdaption(const YAML::Node &adapt_manager)
Initialise adaptive MCMC.
TMatrixDSym * throwMatrix
Matrix which we use for step proposal before Cholesky decomposition (not actually used for step propo...
void ThrowParCurr(const double mag=1.)
Helper function to throw the current parameter by mag sigma. Can study bias in MCMC with this; put di...
double _fGlobalStepScale
Global step scale applied to all params in this class.
void ReserveMemory(const int size)
Initialise vectors with parameters information.
std::vector< std::string > _fFancyNames
Fancy name for example rather than xsec_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.
void ThrowParameters()
Throw the parameters according to the covariance matrix. This shouldn't be used in MCMC code ase it c...
void Init(const std::string &name, const std::string &file)
Initialisation of the class using matrix from root file.
std::string matrixName
Name of cov matrix.
void MakeClosestPosDef(TMatrixDSym *cov)
HW: Finds closest possible positive definite matrix in Frobenius Norm ||.||_frob Where ||X||_frob=sqr...
void RandomConfiguration()
Randomly throw the parameters in their 1 sigma range.
std::vector< double > FlipParameterPoint
Central points around which parameters are flipped.
void UpdateAdaptiveCovariance()
Method to update adaptive MCMC .
std::vector< double > _fError
Prior error on the parameter.
void PrintNominal() const
Print prior value for every parameter.
YAML::Node _fYAMLDoc
Stores config describing systematics.
bool IsParameterFixed(const int i) const
Is parameter fixed or not.
void UpdateThrowMatrix(TMatrixDSym *cov)
void CircularParBounds(const int i, const double LowBound, const double UpBound)
HW :: This method is a tad hacky but modular arithmetic gives me a headache.
void SaveUpdatedMatrixConfig()
KS: After step scale, prefit etc. value were modified save this modified config.
void ResetIndivStepScale()
Adaptive Step Tuning Stuff.
std::unique_ptr< PCAHandler > PCAObj
Struct containing information about PCA.
int _fNumPar
Number of systematic parameters.
double CalcLikelihood() const _noexcept_
Calc penalty term based on inverted covariance matrix.
void FlipParameterValue(const int index, const double FlipPoint)
KS: Flip parameter around given value, for example mass ordering around 0.
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.
void PrintNominalCurrProp() const
Print prior, current and proposed value for each parameter.
ParameterHandlerBase(const std::vector< std::string > &YAMLFile, std::string name, double threshold=-1, int FirstPCAdpar=-999, int LastPCAdpar=-999)
ETA - constructor for a YAML file.
std::vector< std::vector< double > > InvertCovMatrix
KS: Same as above but much faster as TMatrixDSym cache miss.
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...
void PrintIndivStepScale() const
Print step scale for each parameter.
std::vector< std::string > _fNames
ETA _fNames is set automatically in the covariance class to be something like xsec_i,...
std::vector< std::pair< double, double > > CircularBoundsValues
Circular bounds for each parameter (lower, upper)
void EnableSpecialProposal(const YAML::Node &param, const int Index)
Enable special proposal.
double * corr_throw
Result of multiplication of Cholesky matrix and randParams.
void MatchMaCh3OutputBranches(TTree *PosteriorFile, std::vector< double > &BranchValues, std::vector< std::string > &BranchNames)
Matches branches in a TTree to parameters in a systematic handler.
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.
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.
void SpecialStepProposal()
Perform Special Step Proposal.
int PrintLength
KS: This is used when printing parameters, sometimes we have super long parameters name,...
void CorrelateSteps() _noexcept_
Use Cholesky throw matrix for better step proposal.
std::vector< double > _fPropVal
Proposed value of the parameter.
std::vector< std::vector< std::string > > _fSampleNames
Tells to which samples object param should be applied.
int CheckBounds() const _noexcept_
Check if parameters were proposed outside physical boundary.
void ToggleFixAllParameters()
fix parameters at prior values
int GetNumParams() const
Get total number of parameters.
std::string GetParName(const int i) const
Get name of parameter.
std::vector< double > GetProposed() const
Get vector of all proposed parameter values.
virtual double GetLikelihood()
Return CalcLikelihood if some params were thrown out of boundary return LARGE_LOGL
int GetParIndex(const std::string &name) const
Get index based on name.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
std::string GetName() const
Get name of covariance.
double GetParInit(const int i) const
Get prior parameter value.
double GetDiagonalError(const int i) const
Get diagonal error for ith parameter.
void SetName(const std::string &name)
Set matrix name.
void SetStepScale(const double scale, const bool verbose=true)
Set global step scale for covariance object.
void SetCovMatrix(TMatrixDSym *cov)
Set covariance matrix.
void SetFlatPrior(const int i, const bool eL)
Set if parameter should have flat prior or not.
void SetThrowMatrix(TMatrixDSym *cov)
Use new throw matrix, used in adaptive MCMC.
void SetParameters(const std::vector< double > &pars={})
Set parameter values using vector, it has to have same size as covariance class.
void SetBranches(TTree &tree, const bool SaveProposal=false)
set branches for output file
void SetParCurrProp(const int i, const double val)
Set current parameter value.
void SetSingleParameter(const int parNo, const double parVal)
Set value of single param to a given value.
void SetPar(const int i, const double val)
Set all the covariance matrix parameters to a user-defined value.
void SetIndivStepScale(const int ParameterIndex, const double StepScale)
DB Function to set fIndivStepScale from a vector (Can be used from execs and inside covariance constr...
void SetTune(const std::string &TuneName)
KS: Set proposed parameter values vector to be base on tune values, for example set proposed values t...
constexpr static const double _LARGE_LOGL_
Large Likelihood is used it parameter go out of physical boundary, this indicates in MCMC that such s...
Definition: Core.h:73
void MatrixVectorMulti(double *_restrict_ VecMulti, double **_restrict_ matrix, const double *_restrict_ vector, const int n)
KS: Custom function to perform multiplication of matrix and vector with multithreading.
int GetThreadIndex()
thread index inside parallel loop
Definition: Monitor.h:86
double MatrixVectorMultiSingle(double **_restrict_ matrix, const double *_restrict_ vector, const int Length, const int i)
KS: Custom function to perform multiplication of matrix and single element which is thread safe.
void MakeMatrixPosDef(TMatrixDSym *cov)
Makes sure that matrix is positive-definite by adding a small number to on-diagonal elements.
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:48
int GetNThreads()
number of threads which we need for example for TRandom3
Definition: Monitor.cpp:357
std::vector< std::vector< double > > GetCholeskyDecomposedMatrix(const TMatrixDSym &matrix, const std::string &matrixName)
Computes Cholesky decomposition of a symmetric positive definite matrix using custom function which c...
std::string FormatDouble(const double value, const int precision)
Convert double into string for precision, useful for playing with yaml if you don't want to have in c...