MaCh3  2.5.0
Reference Guide
Public Member Functions | Private Member Functions | Private Attributes | List of all members
AdaptiveMCMCHandler Class Reference

Contains information about adaptive covariance matrix [15] [10]. More...

#include <Parameters/AdaptiveMCMCHandler.h>

Collaboration diagram for AdaptiveMCMCHandler:
[legend]

Public Member Functions

 AdaptiveMCMCHandler ()
 Constructor. More...
 
virtual ~AdaptiveMCMCHandler ()
 Destructor. More...
 
void Print () const
 Print all class members. More...
 
bool InitFromConfig (const YAML::Node &adapt_manager, const std::string &matrix_name_str, const std::vector< std::string > *parameter_names, const std::vector< double > *parameters, const std::vector< double > *fixed, const std::vector< bool > *param_skip_adapt_flags, const TMatrixDSym *throwMatrix, const double initial_scale_)
 Read initial values from config file. More...
 
void CreateNewAdaptiveCovariance ()
 If we don't have a covariance matrix to start from for adaptive tune we need to make one! More...
 
void SetAdaptiveBlocks (const std::vector< std::vector< int >> &block_indices)
 HW: sets adaptive block matrix. More...
 
void SaveAdaptiveToFile (const std::string &outFileName, const std::string &systematicName, const bool is_final=false)
 HW: Save adaptive throw matrix to file. More...
 
void SetThrowMatrixFromFile (const std::string &matrix_file_name, const std::string &matrix_name, const std::string &means_name, bool &use_adaptive)
 sets throw matrix from a file More...
 
void CheckMatrixValidityForAdaption (const TMatrixDSym *covMatrix) const
 Check if there are structures in matrix that could result in failure of adaption fits. this mostly cover cases of highly correlated block which seem to results in adaption failure. More...
 
void UpdateAdaptiveCovariance ()
 Method to update adaptive MCMC [15]. More...
 
bool IndivStepScaleAdapt () const
 Tell whether we want reset step scale or not. More...
 
bool UpdateMatrixAdapt ()
 Tell whether matrix should be updated. More...
 
bool AdaptionUpdate () const
 To be fair not a clue... More...
 
bool SkipAdaption () const
 Tell if we are Skipping Adaption. More...
 
int GetParIndex (const std::string &name) const
 Get the index of the parameter given its name. More...
 
std::vector< int > GetParIndicesFromNames (const std::vector< std::string > &param_names) const
 Convert vector of parameter names to indices. More...
 
void SetParams (const std::vector< double > *params)
 Set the current values of the parameters. More...
 
void SetFixed (const std::vector< double > *fix)
 Set the fixed parameters. More...
 
int GetNFixed () const
 Get the number of fixed parameters. More...
 
int GetNumParams () const
 Get the current values of the parameters. More...
 
bool IsFixed (const int ipar) const
 Check if a parameter is fixed. More...
 
double CurrVal (const int par_index) const
 Get Current value of parameter. More...
 
int GetTotalSteps () const
 Get Total Number of Steps. More...
 
void SetTotalSteps (const int nsteps)
 Change Total Number of Steps to new value. More...
 
void IncrementNSteps ()
 Increase by one number of total steps. More...
 
void IncrementAcceptedSteps ()
 
TMatrixDSym * GetAdaptiveCovariance () const
 Increase by one number of total steps. More...
 
std::vector< double > GetParameterMeans () const
 Get the parameter means used in the adaptive handler. More...
 
std::string GetOutFileName () const
 Get Name of Output File. More...
 
double GetAdaptionScale ()
 Get the current adaption scale. More...
 
bool GetUseRobbinsMonro () const
 Use Robbins-Monro approach? More...
 
void UpdateRobbinsMonroScale ()
 Update the scale factor for Robbins-Monro adaption. More...
 

Private Member Functions

void CalculateRobbinsMonroStepLength ()
 Calculate the constant step length for Robbins-Monro adaption. More...
 

Private Attributes

int start_adaptive_throw
 
int start_adaptive_update
 When do we stop update the adaptive matrix. More...
 
int end_adaptive_update
 Steps between changing throw matrix. More...
 
int adaptive_update_step
 Steps between changing throw matrix. More...
 
int adaptive_save_n_iterations
 
std::string output_file_name
 Name of the file to save the adaptive matrices into. More...
 
std::vector< int > adapt_block_matrix_indices
 Indices for block-matrix adaption. More...
 
std::vector< int > adapt_block_sizes
 Size of blocks for adaption. More...
 
std::vector< double > par_means
 Mean values for all parameters. More...
 
TMatrixDSym * adaptive_covariance
 Full adaptive covariance matrix. More...
 
int total_steps
 Total number of MCMC steps. More...
 
double adaption_scale
 Scaling factor. More...
 
double initial_scale
 Initial scaling factor. More...
 
const std::vector< double > * _fFixedPars
 Vector of fixed parameters. More...
 
const std::vector< double > * _fCurrVal
 Current values of parameters. More...
 
const std::vector< std::string > * _ParamNames
 Parameter names. More...
 
int acceptance_rate_batch_size
 Acceptance rate in the current batch. More...
 
bool use_robbins_monro
 Use Robbins Monro https://arxiv.org/pdf/1006.3690. More...
 
double target_acceptance
 Target acceptance rate for Robbins Monro. More...
 
double c_robbins_monro
 Constant "step scaling" factor for Robbins-Monro. More...
 
int n_rm_restarts
 Number of restarts for Robbins Monro (so far) More...
 
int total_rm_restarts
 Total number of restarts ALLOWED for Robbins Monro. More...
 
bool prev_step_accepted
 Need to keep track of whether previous step was accepted for RM. More...
 
const std::vector< bool > * _param_skip_adapt_flags
 Parameters to skip during adaption. More...
 
const TMatrixDSym * initial_throw_matrix
 Initial throw matrix. More...
 
bool initial_throw_matrix_saved = false
 Flag to check if initial throw matrix has been saved. More...
 

Detailed Description

Contains information about adaptive covariance matrix [15] [10].

Author
Henry Wallace
Hank Hua

Adaptation can be performed in two ways:

Definition at line 25 of file AdaptiveMCMCHandler.h.

Constructor & Destructor Documentation

◆ AdaptiveMCMCHandler()

AdaptiveMCMCHandler::AdaptiveMCMCHandler ( )

Constructor.

Robbins-Monro adaption

Definition at line 4 of file AdaptiveMCMCHandler.cpp.

4  {
5 // ********************************************
10  total_steps = 0;
11  par_means = {};
12  adaptive_covariance = nullptr;
13 
15  use_robbins_monro = false;
17 
18  target_acceptance = 0.234;
20 }
int end_adaptive_update
Steps between changing throw matrix.
int total_steps
Total number of MCMC steps.
double target_acceptance
Target acceptance rate for Robbins Monro.
int adaptive_update_step
Steps between changing throw matrix.
int start_adaptive_update
When do we stop update the adaptive matrix.
TMatrixDSym * adaptive_covariance
Full adaptive covariance matrix.
int acceptance_rate_batch_size
Acceptance rate in the current batch.
std::vector< double > par_means
Mean values for all parameters.
int total_rm_restarts
Total number of restarts ALLOWED for Robbins Monro.
bool use_robbins_monro
Use Robbins Monro https://arxiv.org/pdf/1006.3690.

◆ ~AdaptiveMCMCHandler()

AdaptiveMCMCHandler::~AdaptiveMCMCHandler ( )
virtual

Destructor.

Definition at line 23 of file AdaptiveMCMCHandler.cpp.

23  {
24 // ********************************************
25  if(adaptive_covariance != nullptr) {
26  delete adaptive_covariance;
27  }
28 }

Member Function Documentation

◆ AdaptionUpdate()

bool AdaptiveMCMCHandler::AdaptionUpdate ( ) const

To be fair not a clue...

Definition at line 426 of file AdaptiveMCMCHandler.cpp.

426  {
427 // ********************************************
428  if(total_steps <= start_adaptive_throw) return true;
429  else return false;
430 }

◆ CalculateRobbinsMonroStepLength()

void AdaptiveMCMCHandler::CalculateRobbinsMonroStepLength ( )
private

Calculate the constant step length for Robbins-Monro adaption.

Firstly we need to calculate the alpha value, this is in some sense "optimal"

Now we can calculate the scale factor

Definition at line 514 of file AdaptiveMCMCHandler.cpp.

514  {
515 // ********************************************
516  /*
517  Obtains the constant "step scale" for Robbins-Monro adaption
518  for simplicity and because it has the degrees of freedom "baked in"
519  it will be same across ALL blocks.
520  */
521 
523  double alpha = -ROOT::Math::normal_quantile(target_acceptance/2, 1.0);
524 
525  int non_fixed_pars = GetNumParams()-GetNFixed();
526 
527  if(non_fixed_pars == 0){
528  MACH3LOG_ERROR("Number of non_fixed_pars is 0, have you fixed all parameters");
529  MACH3LOG_ERROR("This will cause division by 0 and crash ()");
530  throw MaCh3Exception(__FILE__, __LINE__);
531  }
533  c_robbins_monro = (1- 1/non_fixed_pars)*std::sqrt(TMath::Pi()*2 * std::exp(alpha*alpha));
534  c_robbins_monro += 1/(non_fixed_pars*target_acceptance*(1-target_acceptance));
535 
536  MACH3LOG_INFO("Robbins-Monro scale factor set to {}", c_robbins_monro);
537 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
int GetNumParams() const
Get the current values of the parameters.
double c_robbins_monro
Constant "step scaling" factor for Robbins-Monro.
int GetNFixed() const
Get the number of fixed parameters.
Custom exception class used throughout MaCh3.

◆ CheckMatrixValidityForAdaption()

void AdaptiveMCMCHandler::CheckMatrixValidityForAdaption ( const TMatrixDSym *  covMatrix) const

Check if there are structures in matrix that could result in failure of adaption fits. this mostly cover cases of highly correlated block which seem to results in adaption failure.

Parameters
covMatrixcovariance matrix which we convert into correlation matrix

Definition at line 448 of file AdaptiveMCMCHandler.cpp.

448  {
449 // ********************************************
450  int n = covMatrix->GetNrows();
451  std::vector<std::vector<double>> corrMatrix(n, std::vector<double>(n));
452 
453  // Extract standard deviations from the diagonal
454  std::vector<double> stdDevs(n);
455  for (int i = 0; i < n; ++i) {
456  stdDevs[i] = std::sqrt((*covMatrix)(i, i));
457  }
458 
459  // Compute correlation for each element
460  for (int i = 0; i < n; ++i) {
461  for (int j = 0; j < n; ++j) {
462  double cov = (*covMatrix)(i, j);
463  double corr = cov / (stdDevs[i] * stdDevs[j]);
464  corrMatrix[i][j] = corr;
465  }
466  }
467 
468  // KS: These numbers are completely arbitrary
469  constexpr double NumberOfParametersThreshold = 5;
470  constexpr double CorrelationThreshold = 0.5;
471 
472  bool FlagMessage = false;
473  // For each parameter, count how many others it is correlated with
474  for (int i = 0; i < n; ++i) {
475  if (IsFixed(i)) {
476  continue;
477  }
478  int block_i = adapt_block_matrix_indices[i];
479  int correlatedCount = 0;
480  std::vector<int> correlatedWith;
481 
482  for (int j = 0; j < n; ++j) {
483  if (i == j || IsFixed(j) || adapt_block_matrix_indices[j] != block_i) {
484  continue;
485  }
486  if (corrMatrix[i][j] > CorrelationThreshold) {
487  correlatedCount++;
488  correlatedWith.push_back(j);
489  }
490  }
491 
492  if (correlatedCount > NumberOfParametersThreshold) {
493  MACH3LOG_WARN("Parameter {} ({}) is highly correlated with {} other parameters (above threshold {}).", _ParamNames->at(i), i, correlatedCount, CorrelationThreshold);
494  MACH3LOG_WARN("Correlated with parameters: {}.", fmt::join(correlatedWith, ", "));
495  FlagMessage = true;
496  }
497  }
498  if(FlagMessage){
499  MACH3LOG_WARN("Found highly correlated block, adaptive may not work as intended, proceed with caution");
500  MACH3LOG_WARN("You can consider defying so called Adaption Block to not update this part of matrix");
501  }
502 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
const std::vector< std::string > * _ParamNames
Parameter names.
bool IsFixed(const int ipar) const
Check if a parameter is fixed.
std::vector< int > adapt_block_matrix_indices
Indices for block-matrix adaption.

◆ CreateNewAdaptiveCovariance()

void AdaptiveMCMCHandler::CreateNewAdaptiveCovariance ( )

If we don't have a covariance matrix to start from for adaptive tune we need to make one!

Definition at line 142 of file AdaptiveMCMCHandler.cpp.

142  {
143 // ********************************************
144  adaptive_covariance = new TMatrixDSym(GetNumParams());
145  adaptive_covariance->Zero();
146  par_means = std::vector<double>(GetNumParams(), 0);
147 }

◆ CurrVal()

double AdaptiveMCMCHandler::CurrVal ( const int  par_index) const

Get Current value of parameter.

HW Implemented as its own method to allow for different behaviour in the future

Definition at line 506 of file AdaptiveMCMCHandler.cpp.

506  {
507 // ********************************************
510  return (*_fCurrVal)[par_index];
511 }
const std::vector< double > * _fCurrVal
Current values of parameters.

◆ GetAdaptionScale()

double AdaptiveMCMCHandler::GetAdaptionScale ( )
inline

Get the current adaption scale.

Definition at line 222 of file AdaptiveMCMCHandler.h.

223  {
224  return adaption_scale;
225  }
double adaption_scale
Scaling factor.

◆ GetAdaptiveCovariance()

TMatrixDSym* AdaptiveMCMCHandler::GetAdaptiveCovariance ( ) const
inline

Increase by one number of total steps.

Definition at line 202 of file AdaptiveMCMCHandler.h.

203  {
204  return adaptive_covariance;
205  }

◆ GetNFixed()

int AdaptiveMCMCHandler::GetNFixed ( ) const
inline

Get the number of fixed parameters.

Definition at line 139 of file AdaptiveMCMCHandler.h.

140  {
141  if (!_fFixedPars)
142  {
143  return 0;
144  }
145  int n_fixed = 0;
146  for (int i = 0; i < static_cast<int>(_fFixedPars->size()); i++)
147  {
148  if ((*_fFixedPars)[i] < 0)
149  {
150  n_fixed++;
151  }
152  }
153  return n_fixed;
154  }
const std::vector< double > * _fFixedPars
Vector of fixed parameters.

◆ GetNumParams()

int AdaptiveMCMCHandler::GetNumParams ( ) const
inline

Get the current values of the parameters.

Definition at line 158 of file AdaptiveMCMCHandler.h.

159  {
160  return static_cast<int>(_fCurrVal->size());
161  }

◆ GetOutFileName()

std::string AdaptiveMCMCHandler::GetOutFileName ( ) const
inline

Get Name of Output File.

Definition at line 216 of file AdaptiveMCMCHandler.h.

217  {
218  return output_file_name;
219  }
std::string output_file_name
Name of the file to save the adaptive matrices into.

◆ GetParameterMeans()

std::vector<double> AdaptiveMCMCHandler::GetParameterMeans ( ) const
inline

Get the parameter means used in the adaptive handler.

Definition at line 209 of file AdaptiveMCMCHandler.h.

210  {
211  return par_means;
212  }

◆ GetParIndex()

int AdaptiveMCMCHandler::GetParIndex ( const std::string &  name) const
inline

Get the index of the parameter given its name.

Parameters
nameName of the parameter
Returns
Index of the parameter

Definition at line 96 of file AdaptiveMCMCHandler.h.

96  {
97  // HH: Adapted from ParameterHandlerBase.cpp
98  int index = M3::_BAD_INT_;
99  for (size_t i = 0; i < _ParamNames->size(); i++) {
100  if(name == _ParamNames->at(i)) {
101  index = static_cast<int>(i);
102  break;
103  }
104  }
105  return index;
106  }
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55

◆ GetParIndicesFromNames()

std::vector<int> AdaptiveMCMCHandler::GetParIndicesFromNames ( const std::vector< std::string > &  param_names) const
inline

Convert vector of parameter names to indices.

Parameters
param_namesVector of parameter names
Returns
Vector of parameter indices

Definition at line 111 of file AdaptiveMCMCHandler.h.

111  {
112  std::vector<int> indices;
113  for (const auto& name : param_names) {
114  int index = GetParIndex(name);
115  if (index != M3::_BAD_INT_) {
116  indices.push_back(index);
117  } else {
118  MACH3LOG_ERROR("Parameter name {} not found in parameter list", name);
119  throw MaCh3Exception(__FILE__, __LINE__);
120  }
121  }
122  return indices;
123  }
int GetParIndex(const std::string &name) const
Get the index of the parameter given its name.

◆ GetTotalSteps()

int AdaptiveMCMCHandler::GetTotalSteps ( ) const
inline

Get Total Number of Steps.

Definition at line 178 of file AdaptiveMCMCHandler.h.

179  {
180  return total_steps;
181  }

◆ GetUseRobbinsMonro()

bool AdaptiveMCMCHandler::GetUseRobbinsMonro ( ) const
inline

Use Robbins-Monro approach?

Definition at line 228 of file AdaptiveMCMCHandler.h.

229  {
230  return use_robbins_monro;
231  }

◆ IncrementAcceptedSteps()

void AdaptiveMCMCHandler::IncrementAcceptedSteps ( )
inline

Definition at line 195 of file AdaptiveMCMCHandler.h.

196  {
197  prev_step_accepted = true;
198  }
bool prev_step_accepted
Need to keep track of whether previous step was accepted for RM.

◆ IncrementNSteps()

void AdaptiveMCMCHandler::IncrementNSteps ( )
inline

Increase by one number of total steps.

Definition at line 190 of file AdaptiveMCMCHandler.h.

191  {
192  total_steps++;
193  }

◆ IndivStepScaleAdapt()

bool AdaptiveMCMCHandler::IndivStepScaleAdapt ( ) const

Tell whether we want reset step scale or not.

Definition at line 399 of file AdaptiveMCMCHandler.cpp.

399  {
400 // ********************************************
401  if(total_steps == start_adaptive_throw) return true;
402  else return false;
403 }

◆ InitFromConfig()

bool AdaptiveMCMCHandler::InitFromConfig ( const YAML::Node &  adapt_manager,
const std::string &  matrix_name_str,
const std::vector< std::string > *  parameter_names,
const std::vector< double > *  parameters,
const std::vector< double > *  fixed,
const std::vector< bool > *  param_skip_adapt_flags,
const TMatrixDSym *  throwMatrix,
const double  initial_scale_ 
)

Read initial values from config file.

Parameters
adapt_managerYAML node containing the configuration (AdaptionOptions).
matrix_name_strName of the covariance matrix block to configure.
parametersPointer to a vector of parameter values (nominal values).
fixedPointer to a vector of fixed parameter values.
Returns
True if adaptive MCMC configuration was successfully initialized, false otherwise.

HW: This is technically wrong, should be across all systematics but will be addressed in a later PR

Definition at line 31 of file AdaptiveMCMCHandler.cpp.

37  {
38 // ********************************************
39  /*
40  * HW: Idea is that adaption can simply read the YAML config
41  * Options :
42  * External Info:
43  * UseExternalMatrix [bool] : Use an external matrix
44  * ExternalMatrixFileName [str] : Name of file containing external info
45  * ExternalMatrixName [str] : Name of external Matrix
46  * ExternalMeansName [str] : Name of external means vector [for updates]
47  *
48  * General Info:
49  * DoAdaption [bool] : Do we want to do adaption?
50  * AdaptionStartThrow [int] : Step we start throwing adaptive matrix from
51  * AdaptionEndUpdate [int] : Step we stop updating adaptive matrix
52  * AdaptionStartUpdate [int] : Do we skip the first N steps?
53  * AdaptionUpdateStep [int] : Number of steps between matrix updates
54  * AdaptionSaveNIterations [int]: You don't have to save every adaptive stage so decide how often you want to save
55  * Adaption blocks [vector<vector<int>>] : Splits the throw matrix into several block matrices
56  * OuputFileName [std::string] : Name of the file that the adaptive matrices will be saved into
57  */
58 
59  // setAdaptionDefaults();
60  if(!adapt_manager["AdaptionOptions"]["Covariance"][matrix_name_str]) {
61  MACH3LOG_WARN("Adaptive Settings not found for {}, this is fine if you don't want adaptive MCMC", matrix_name_str);
62  return false;
63  }
64 
65  // We"re going to grab this info from the YAML manager
66  if(!GetFromManager<bool>(adapt_manager["AdaptionOptions"]["Covariance"][matrix_name_str]["DoAdaption"], false)) {
67  MACH3LOG_WARN("Not using adaption for {}", matrix_name_str);
68  return false;
69  }
70 
71  if(!CheckNodeExists(adapt_manager, "AdaptionOptions", "Settings", "OutputFileName")) {
72  MACH3LOG_ERROR("No OutputFileName specified in AdaptionOptions::Settings into your config file");
73  MACH3LOG_ERROR("This is required if you are using adaptive MCMC");
74  throw MaCh3Exception(__FILE__, __LINE__);
75  }
76 
77  start_adaptive_throw = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["StartThrow"], 10);
78  start_adaptive_update = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["StartUpdate"], 0);
79  end_adaptive_update = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["EndUpdate"], 10000);
80  adaptive_update_step = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["UpdateStep"], 100);
81  adaptive_save_n_iterations = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["SaveNIterations"], -1);
82  output_file_name = GetFromManager<std::string>(adapt_manager["AdaptionOptions"]["Settings"]["OutputFileName"], "");
83 
84  // Check for Robbins-Monro adaption
85  use_robbins_monro = GetFromManager<bool>(adapt_manager["AdaptionOptions"]["Settings"]["UseRobbinsMonro"], false);
86  target_acceptance = GetFromManager<double>(adapt_manager["AdaptionOptions"]["Settings"]["TargetAcceptance"], 0.234);
87  if (target_acceptance <= 0 || target_acceptance >= 1) {
88  MACH3LOG_ERROR("Target acceptance must be in (0,1), got {}", target_acceptance);
89  throw MaCh3Exception(__FILE__, __LINE__);
90  }
91 
92  acceptance_rate_batch_size = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["AcceptanceRateBatchSize"], 10000);
93  total_rm_restarts = GetFromManager<int>(adapt_manager["AdaptionOptions"]["Settings"]["TotalRobbinsMonroRestarts"], 0);
94  n_rm_restarts = 0;
95 
96  prev_step_accepted = false;
97 
98  // We also want to check for "blocks" by default all parameters "know" about each other
99  // but we can split the matrix into independent block matrices
100  SetParams(parameters);
101  SetFixed(fixed);
102 
103  _ParamNames = parameter_names;
104 
105  initial_throw_matrix = static_cast<TMatrixDSym*>(throwMatrix->Clone());
106 
107  initial_scale = initial_scale_;
108 
110  if(use_robbins_monro) {
111  MACH3LOG_INFO("Using Robbins-Monro for adaptive step size tuning with target acceptance {}", target_acceptance);
112  // Now we need some thing
114  }
115 
116  // We'll use the same start scale for Robbins-Monro as standard adaption
117  adaption_scale = 2.38/std::sqrt(GetNumParams());
118 
119  // HH: added ability to define blocks by parameter names
120  bool matrix_blocks_by_name = GetFromManager<bool>(adapt_manager["AdaptionOptions"]["Covariance"][matrix_name_str]["MatrixBlocksByName"], false);
121  std::vector<std::vector<int>> matrix_blocks;
122  // HH: read blocks by names and convert to indices
123  if (matrix_blocks_by_name) {
124  auto matrix_blocks_input = GetFromManager<std::vector<std::vector<std::string>>>(adapt_manager["AdaptionOptions"]["Covariance"][matrix_name_str]["MatrixBlocks"], {{}});
125  // Now we need to convert to indices
126  for (const auto& block : matrix_blocks_input) {
127  auto block_indices = GetParIndicesFromNames(block);
128  matrix_blocks.push_back(block_indices);
129  }
130  } else {
131  matrix_blocks = GetFromManager<std::vector<std::vector<int>>>(adapt_manager["AdaptionOptions"]["Covariance"][matrix_name_str]["MatrixBlocks"], {{}});
132  }
133  SetAdaptiveBlocks(matrix_blocks);
134 
135  // set parameters to skip during adaption
136  _param_skip_adapt_flags = param_skip_adapt_flags;
137 
138  return true;
139 }
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60
const std::vector< bool > * _param_skip_adapt_flags
Parameters to skip during adaption.
std::vector< int > GetParIndicesFromNames(const std::vector< std::string > &param_names) const
Convert vector of parameter names to indices.
void SetFixed(const std::vector< double > *fix)
Set the fixed parameters.
const TMatrixDSym * initial_throw_matrix
Initial throw matrix.
void SetParams(const std::vector< double > *params)
Set the current values of the parameters.
double initial_scale
Initial scaling factor.
void CalculateRobbinsMonroStepLength()
Calculate the constant step length for Robbins-Monro adaption.
void SetAdaptiveBlocks(const std::vector< std::vector< int >> &block_indices)
HW: sets adaptive block matrix.
int n_rm_restarts
Number of restarts for Robbins Monro (so far)

◆ IsFixed()

bool AdaptiveMCMCHandler::IsFixed ( const int  ipar) const
inline

Check if a parameter is fixed.

Definition at line 164 of file AdaptiveMCMCHandler.h.

165  {
166  if (!_fFixedPars)
167  {
168  return false;
169  }
170  return ((*_fFixedPars)[ipar] < 0);
171  }

◆ Print()

void AdaptiveMCMCHandler::Print ( ) const

Print all class members.

Definition at line 433 of file AdaptiveMCMCHandler.cpp.

433  {
434 // ********************************************
435  MACH3LOG_INFO("Adaptive MCMC Info:");
436  MACH3LOG_INFO("Throwing from New Matrix from Step : {}", start_adaptive_throw);
437  MACH3LOG_INFO("Adaption Matrix Start Update : {}", start_adaptive_update);
438  MACH3LOG_INFO("Adaption Matrix Ending Updates : {}", end_adaptive_update);
439  MACH3LOG_INFO("Steps Between Updates : {}", adaptive_update_step);
440  MACH3LOG_INFO("Saving matrices to file : {}", output_file_name);
441  MACH3LOG_INFO("Will only save every {} iterations" , adaptive_save_n_iterations);
442  if(use_robbins_monro) {
443  MACH3LOG_INFO("Using Robbins-Monro for adaptive step size tuning with target acceptance {}", target_acceptance);
444  }
445 }

◆ SaveAdaptiveToFile()

void AdaptiveMCMCHandler::SaveAdaptiveToFile ( const std::string &  outFileName,
const std::string &  systematicName,
const bool  is_final = false 
)

HW: Save adaptive throw matrix to file.

Definition at line 192 of file AdaptiveMCMCHandler.cpp.

193  {
194 // ********************************************
195  // Skip saving if adaptive_save_n_iterations is negative,
196  // unless this is the final iteration (is_final overrides the condition)
197  if (adaptive_save_n_iterations < 0 && !is_final) return;
198  if (is_final ||
201  TFile *outFile = new TFile(outFileName.c_str(), "UPDATE");
202  if (outFile->IsZombie()) {
203  MACH3LOG_ERROR("Couldn't find {}", outFileName);
204  throw MaCh3Exception(__FILE__, __LINE__);
205  }
206 
207  TVectorD *outMeanVec = new TVectorD(int(par_means.size()));
208  for (int i = 0; i < int(par_means.size()); i++) {
209  (*outMeanVec)(i) = par_means[i];
210  }
211 
212  std::string adaptive_cov_name = systematicName + "_posfit_matrix";
213  std::string mean_vec_name = systematicName + "_mean_vec";
214  if (!is_final) {
215  // Some string to make the name of the saved adaptive matrix clear
216  std::string total_steps_str =
217  std::to_string(total_steps);
218  std::string syst_name_str = systematicName;
219  adaptive_cov_name =
220  total_steps_str + '_' + syst_name_str + std::string("_throw_matrix");
221  mean_vec_name =
222  total_steps_str + '_' + syst_name_str + std::string("_mean_vec");
223  }
224 
225  outFile->cd();
226  auto adaptive_to_save = static_cast<TMatrixDSym*>(adaptive_covariance->Clone());
227  // Multiply by scale
228  // HH: Only scale parameters not being skipped
229  for (int i = 0; i < GetNumParams(); ++i) {
230  for (int j = 0; j <= i; ++j) {
231  // If both parameters are skipped, scale by initial scale
232  if ((*_param_skip_adapt_flags)[i] && (*_param_skip_adapt_flags)[j]) {
233  (*adaptive_to_save)(i, j) *= initial_scale * initial_scale;
234  if (i != j) {
235  (*adaptive_to_save)(j, i) *= initial_scale * initial_scale;
236  }
237  }
238  else if (!(*_param_skip_adapt_flags)[i] && !(*_param_skip_adapt_flags)[j]) {
239  (*adaptive_to_save)(i, j) *= GetAdaptionScale() * GetAdaptionScale();
240  if (i != j) {
241  (*adaptive_to_save)(j, i) *= GetAdaptionScale() * GetAdaptionScale();
242  }
243  }
244  else {
245  // Not too sure what to do here, as the correlation should be zero
246  // Scale by both factors as a compromise
247  (*adaptive_to_save)(i, j) *= GetAdaptionScale() * initial_scale;
248  if (i != j) {
249  (*adaptive_to_save)(j, i) *= GetAdaptionScale() * initial_scale;
250  }
251  }
252  }
253  }
254 
255  adaptive_to_save->Write(adaptive_cov_name.c_str());
256  outMeanVec->Write(mean_vec_name.c_str());
257 
258  // HH: Also save the initial throw matrix for reference
260  std::string initial_cov_name = "0_" + systematicName + std::string("_initial_throw_matrix");
261  auto initial_throw_matrix_to_save = static_cast<TMatrixDSym*>(initial_throw_matrix->Clone());
262  *(initial_throw_matrix_to_save) *= initial_scale * initial_scale;
264  initial_throw_matrix_to_save->Write(initial_cov_name.c_str());
265  delete initial_throw_matrix_to_save;
266  }
267 
268  outFile->Close();
269  delete adaptive_to_save;
270  delete outMeanVec;
271  delete outFile;
272  }
273 }
double GetAdaptionScale()
Get the current adaption scale.
bool initial_throw_matrix_saved
Flag to check if initial throw matrix has been saved.

◆ SetAdaptiveBlocks()

void AdaptiveMCMCHandler::SetAdaptiveBlocks ( const std::vector< std::vector< int >> &  block_indices)

HW: sets adaptive block matrix.

Parameters
block_indicesValues for sub-matrix blocks

Definition at line 150 of file AdaptiveMCMCHandler.cpp.

150  {
151 // ********************************************
152  /*
153  * In order to adapt efficient we want to setup our throw matrix to be a serious of block-diagonal (ish) matrices
154  *
155  * To do this we set sub-block in the config by parameter index. For example having
156  * [[0,4],[4, 6]] in your config will set up two blocks one with all indices 0<=i<4 and the other with 4<=i<6
157  */
158  // Set up block regions
159  adapt_block_matrix_indices = std::vector<int>(GetNumParams(), 0);
160 
161  // Should also make a matrix of block sizes
162  adapt_block_sizes = std::vector<int>(block_indices.size()+1, 0);
164 
165  if(block_indices.size()==0 || block_indices[0].size()==0) return;
166 
167  int block_size = static_cast<int>(block_indices.size());
168  // Now we loop over our blocks
169  for(int iblock=0; iblock < block_size; iblock++){
170  // Loop over blocks in the block
171  int sub_block_size = static_cast<int>(block_indices.size()-1);
172  for(int isubblock=0; isubblock < sub_block_size ; isubblock+=2){
173  int block_lb = block_indices[iblock][isubblock];
174  int block_ub = block_indices[iblock][isubblock+1];
175 
176  if(block_lb > GetNumParams() || block_ub > GetNumParams()){
177  MACH3LOG_ERROR("Cannot set matrix block with edges {}, {} for matrix of size {}",
178  block_lb, block_ub, GetNumParams());
179  throw MaCh3Exception(__FILE__, __LINE__);;
180  }
181  for(int ipar = block_lb; ipar < block_ub; ipar++){
182  adapt_block_matrix_indices[ipar] = iblock+1;
183  adapt_block_sizes[iblock+1] += 1;
184  adapt_block_sizes[0] -= 1;
185  }
186  }
187  }
188 }
std::vector< int > adapt_block_sizes
Size of blocks for adaption.

◆ SetFixed()

void AdaptiveMCMCHandler::SetFixed ( const std::vector< double > *  fix)
inline

Set the fixed parameters.

Definition at line 132 of file AdaptiveMCMCHandler.h.

133  {
134  _fFixedPars = fix;
135  }

◆ SetParams()

void AdaptiveMCMCHandler::SetParams ( const std::vector< double > *  params)
inline

Set the current values of the parameters.

Definition at line 126 of file AdaptiveMCMCHandler.h.

127  {
128  _fCurrVal = params;
129  }

◆ SetThrowMatrixFromFile()

void AdaptiveMCMCHandler::SetThrowMatrixFromFile ( const std::string &  matrix_file_name,
const std::string &  matrix_name,
const std::string &  means_name,
bool &  use_adaptive 
)

sets throw matrix from a file

Parameters
matrix_file_namename of file matrix lives in
matrix_namename of matrix in file
means_namename of means vec in file

Definition at line 278 of file AdaptiveMCMCHandler.cpp.

281  {
282 // ********************************************
283  // Lets you set the throw matrix externally
284  // Open file
285  auto matrix_file = std::make_unique<TFile>(matrix_file_name.c_str());
286  use_adaptive = true;
287 
288  if(matrix_file->IsZombie()){
289  MACH3LOG_ERROR("Couldn't find {}", matrix_file_name);
290  throw MaCh3Exception(__FILE__ , __LINE__ );
291  }
292 
293  // Next we grab our matrix
294  adaptive_covariance = static_cast<TMatrixDSym*>(matrix_file->Get(matrix_name.c_str()));
295  if(!adaptive_covariance){
296  MACH3LOG_ERROR("Couldn't find {} in {}", matrix_name, matrix_file_name);
297  throw MaCh3Exception(__FILE__ , __LINE__ );
298  }
299 
300  // Finally we grab the means vector
301  TVectorD* means_vector = static_cast<TVectorD*>(matrix_file->Get(means_name.c_str()));
302 
303  // This is fine to not exist!
304  if(means_vector){
305  // Yay our vector exists! Let's loop and fill it
306  // Should check this is done
307  if(means_vector->GetNrows() != GetNumParams()){
308  MACH3LOG_ERROR("External means vec size ({}) != matrix size ({})", means_vector->GetNrows(), GetNumParams());
309  throw MaCh3Exception(__FILE__, __LINE__);
310  }
311 
312  par_means = std::vector<double>(GetNumParams());
313  for(int i = 0; i < GetNumParams(); i++){
314  par_means[i] = (*means_vector)(i);
315  }
316  MACH3LOG_INFO("Found Means in External File, Will be able to adapt");
317  }
318  // Totally fine if it doesn't exist, we just can't do adaption
319  else{
320  // We don't need a means vector, set the adaption=false
321  MACH3LOG_WARN("Cannot find means vector in {}, therefore I will not be able to adapt!", matrix_file_name);
322  use_adaptive = false;
323  }
324 
325  matrix_file->Close();
326  MACH3LOG_INFO("Set up matrix from external file");
327 }

◆ SetTotalSteps()

void AdaptiveMCMCHandler::SetTotalSteps ( const int  nsteps)
inline

Change Total Number of Steps to new value.

Definition at line 184 of file AdaptiveMCMCHandler.h.

185  {
186  total_steps = nsteps;
187  }

◆ SkipAdaption()

bool AdaptiveMCMCHandler::SkipAdaption ( ) const

Tell if we are Skipping Adaption.

Definition at line 418 of file AdaptiveMCMCHandler.cpp.

418  {
419 // ********************************************
421  total_steps< start_adaptive_update) return true;
422  else return false;
423 }

◆ UpdateAdaptiveCovariance()

void AdaptiveMCMCHandler::UpdateAdaptiveCovariance ( )

Method to update adaptive MCMC [15].

Parameters
_fCurrValValue of each parameter necessary for updating throw matrix

Definition at line 330 of file AdaptiveMCMCHandler.cpp.

330  {
331 // ********************************************
332  std::vector<double> par_means_prev = par_means;
333  int steps_post_burn = total_steps - start_adaptive_update;
334 
335  prev_step_accepted = false;
336 
337  // Step 1: Update means and compute deviations
338  for (int i = 0; i < GetNumParams(); ++i) {
339  if (IsFixed(i)) continue;
340 
341  par_means[i] = (CurrVal(i) + par_means_prev[i]*steps_post_burn)/(steps_post_burn+1);
342  // Left over from cyclic means
343  }
344 
345  // Step 2: Update covariance
346  #ifdef MULTITHREAD
347  #pragma omp parallel for
348  #endif
349  for (int i = 0; i < GetNumParams(); ++i) {
350  if (IsFixed(i)) {
351  (*adaptive_covariance)(i, i) = 1.0;
352  continue;
353  }
354 
355  int block_i = adapt_block_matrix_indices[i];
356 
357  for (int j = 0; j <= i; ++j) {
358  if (IsFixed(j) || adapt_block_matrix_indices[j] != block_i) {
359  (*adaptive_covariance)(i, j) = 0.0;
360  (*adaptive_covariance)(j, i) = 0.0;
361  continue;
362  }
363 
364  double cov_prev = (*adaptive_covariance)(i, j);
365  double cov_updated = 0.0;
366 
367  // HH: Check if either parameter is skipped
368  // If parameter i and j are both skipped
369  if ((*_param_skip_adapt_flags)[i] && (*_param_skip_adapt_flags)[j]) {
370  // Set covariance to initial throw matrix value
371  (*adaptive_covariance)(i, j) = (*initial_throw_matrix)(i, j);
372  (*adaptive_covariance)(j, i) = (*initial_throw_matrix)(j, i);
373  continue;
374  }
375  // If only one parameter is skipped we need to make sure correlation is zero
376  // because we don't want to change one parameter while keeping the other fixed
377  // as this would lead to penalty terms in the prior
378  else if ((*_param_skip_adapt_flags)[i] || (*_param_skip_adapt_flags)[j]) {
379  continue;
380  } // otherwise adapt as normal
381 
382  if (steps_post_burn > 0) {
383  // Haario-style update
384  double cov_t = cov_prev * (steps_post_burn - 1) / steps_post_burn;
385  double prev_means_t = steps_post_burn * par_means_prev[i] * par_means_prev[j];
386  double curr_means_t = (steps_post_burn + 1) * par_means[i] * par_means[j];
387  double curr_step_t = CurrVal(i) * CurrVal(j);
388 
389  cov_updated = cov_t + (prev_means_t - curr_means_t + curr_step_t) / steps_post_burn;
390  }
391 
392  (*adaptive_covariance)(i, j) = cov_updated;
393  (*adaptive_covariance)(j, i) = cov_updated;
394  }
395  }
396 }
double CurrVal(const int par_index) const
Get Current value of parameter.

◆ UpdateMatrixAdapt()

bool AdaptiveMCMCHandler::UpdateMatrixAdapt ( )

Tell whether matrix should be updated.

Definition at line 406 of file AdaptiveMCMCHandler.cpp.

406  {
407 // ********************************************
409  // Check whether the number of steps is divisible by the adaptive update step
410  // e.g. if adaptive_update_step = 1000 and (total_step - start_adpative_throw) is 5000 then this is true
412  return true;
413  }
414  else return false;
415 }

◆ UpdateRobbinsMonroScale()

void AdaptiveMCMCHandler::UpdateRobbinsMonroScale ( )

Update the scale factor for Robbins-Monro adaption.

This allows to move towards a STABLE number of steps in the batch

Update the scale factor for Robbins-Monro adaption [200 is arbitrary for early stability but motivated by paper]

Now we either increase or decrease the scale

Definition at line 540 of file AdaptiveMCMCHandler.cpp.

540  {
541 // ********************************************
542  /*
543  Update the scale factor using Robbins-Monro.
544  TLDR: If acceptance rate is too high, scale factor goes up, if too low goes down
545  will pretty rapidly converge to the right value.
546  */
547  if ((total_steps>0) && (total_steps % acceptance_rate_batch_size == 0)) {
548  n_rm_restarts++;
549  }
550 
552  int batch_step;
554  batch_step = total_steps;
555  }
556  else{
558  }
559 
560  int non_fixed_pars = GetNumParams()-GetNFixed();
561 
562  // double root_scale = std::sqrt(adaption_scale);
563 
565  double scale_factor = adaption_scale * c_robbins_monro / std::max(200.0, static_cast<double>(batch_step) / non_fixed_pars);
566 
568  if(prev_step_accepted){
569  adaption_scale += (1 - target_acceptance) * scale_factor;
570  }
571  else{
572  adaption_scale -= target_acceptance * scale_factor;
573  }
574 }

Member Data Documentation

◆ _fCurrVal

const std::vector<double>* AdaptiveMCMCHandler::_fCurrVal
private

Current values of parameters.

Definition at line 286 of file AdaptiveMCMCHandler.h.

◆ _fFixedPars

const std::vector<double>* AdaptiveMCMCHandler::_fFixedPars
private

Vector of fixed parameters.

Definition at line 283 of file AdaptiveMCMCHandler.h.

◆ _param_skip_adapt_flags

const std::vector<bool>* AdaptiveMCMCHandler::_param_skip_adapt_flags
private

Parameters to skip during adaption.

Definition at line 313 of file AdaptiveMCMCHandler.h.

◆ _ParamNames

const std::vector<std::string>* AdaptiveMCMCHandler::_ParamNames
private

Parameter names.

Definition at line 289 of file AdaptiveMCMCHandler.h.

◆ acceptance_rate_batch_size

int AdaptiveMCMCHandler::acceptance_rate_batch_size
private

Acceptance rate in the current batch.

Definition at line 292 of file AdaptiveMCMCHandler.h.

◆ adapt_block_matrix_indices

std::vector<int> AdaptiveMCMCHandler::adapt_block_matrix_indices
private

Indices for block-matrix adaption.

Definition at line 261 of file AdaptiveMCMCHandler.h.

◆ adapt_block_sizes

std::vector<int> AdaptiveMCMCHandler::adapt_block_sizes
private

Size of blocks for adaption.

Definition at line 264 of file AdaptiveMCMCHandler.h.

◆ adaption_scale

double AdaptiveMCMCHandler::adaption_scale
private

Scaling factor.

Definition at line 277 of file AdaptiveMCMCHandler.h.

◆ adaptive_covariance

TMatrixDSym* AdaptiveMCMCHandler::adaptive_covariance
private

Full adaptive covariance matrix.

Definition at line 271 of file AdaptiveMCMCHandler.h.

◆ adaptive_save_n_iterations

int AdaptiveMCMCHandler::adaptive_save_n_iterations
private

If you don't want to save every adaption then you can specify this here

Definition at line 255 of file AdaptiveMCMCHandler.h.

◆ adaptive_update_step

int AdaptiveMCMCHandler::adaptive_update_step
private

Steps between changing throw matrix.

Definition at line 251 of file AdaptiveMCMCHandler.h.

◆ c_robbins_monro

double AdaptiveMCMCHandler::c_robbins_monro
private

Constant "step scaling" factor for Robbins-Monro.

Definition at line 301 of file AdaptiveMCMCHandler.h.

◆ end_adaptive_update

int AdaptiveMCMCHandler::end_adaptive_update
private

Steps between changing throw matrix.

Definition at line 248 of file AdaptiveMCMCHandler.h.

◆ initial_scale

double AdaptiveMCMCHandler::initial_scale
private

Initial scaling factor.

Definition at line 280 of file AdaptiveMCMCHandler.h.

◆ initial_throw_matrix

const TMatrixDSym* AdaptiveMCMCHandler::initial_throw_matrix
private

Initial throw matrix.

Definition at line 316 of file AdaptiveMCMCHandler.h.

◆ initial_throw_matrix_saved

bool AdaptiveMCMCHandler::initial_throw_matrix_saved = false
private

Flag to check if initial throw matrix has been saved.

Definition at line 319 of file AdaptiveMCMCHandler.h.

◆ n_rm_restarts

int AdaptiveMCMCHandler::n_rm_restarts
private

Number of restarts for Robbins Monro (so far)

Definition at line 304 of file AdaptiveMCMCHandler.h.

◆ output_file_name

std::string AdaptiveMCMCHandler::output_file_name
private

Name of the file to save the adaptive matrices into.

Definition at line 258 of file AdaptiveMCMCHandler.h.

◆ par_means

std::vector<double> AdaptiveMCMCHandler::par_means
private

Mean values for all parameters.

Definition at line 268 of file AdaptiveMCMCHandler.h.

◆ prev_step_accepted

bool AdaptiveMCMCHandler::prev_step_accepted
private

Need to keep track of whether previous step was accepted for RM.

Definition at line 310 of file AdaptiveMCMCHandler.h.

◆ start_adaptive_throw

int AdaptiveMCMCHandler::start_adaptive_throw
private

Meta variables related to adaption run time When do we start throwing

Definition at line 242 of file AdaptiveMCMCHandler.h.

◆ start_adaptive_update

int AdaptiveMCMCHandler::start_adaptive_update
private

When do we stop update the adaptive matrix.

Definition at line 245 of file AdaptiveMCMCHandler.h.

◆ target_acceptance

double AdaptiveMCMCHandler::target_acceptance
private

Target acceptance rate for Robbins Monro.

Definition at line 298 of file AdaptiveMCMCHandler.h.

◆ total_rm_restarts

int AdaptiveMCMCHandler::total_rm_restarts
private

Total number of restarts ALLOWED for Robbins Monro.

Definition at line 307 of file AdaptiveMCMCHandler.h.

◆ total_steps

int AdaptiveMCMCHandler::total_steps
private

Total number of MCMC steps.

Definition at line 274 of file AdaptiveMCMCHandler.h.

◆ use_robbins_monro

bool AdaptiveMCMCHandler::use_robbins_monro
private

Use Robbins Monro https://arxiv.org/pdf/1006.3690.

Definition at line 295 of file AdaptiveMCMCHandler.h.


The documentation for this class was generated from the following files: