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