MaCh3  2.4.2
Reference Guide
Public Member Functions | Private Member Functions | Private Attributes | List of all members
adaptive_mcmc::AdaptiveMCMCHandler Class Reference

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

#include <Parameters/AdaptiveMCMCHandler.h>

Collaboration diagram for adaptive_mcmc::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 [13]. 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 [13] [10].

Author
Henry Wallace
Hank Hua

Adaptation can be performed in two ways:

Definition at line 26 of file AdaptiveMCMCHandler.h.

Constructor & Destructor Documentation

◆ AdaptiveMCMCHandler()

adaptive_mcmc::AdaptiveMCMCHandler::AdaptiveMCMCHandler ( )

Constructor.

Robbins-Monro adaption

Definition at line 6 of file AdaptiveMCMCHandler.cpp.

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

◆ ~AdaptiveMCMCHandler()

adaptive_mcmc::AdaptiveMCMCHandler::~AdaptiveMCMCHandler ( )
virtual

Destructor.

Definition at line 25 of file AdaptiveMCMCHandler.cpp.

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

Member Function Documentation

◆ AdaptionUpdate()

bool adaptive_mcmc::AdaptiveMCMCHandler::AdaptionUpdate ( ) const

To be fair not a clue...

Definition at line 428 of file AdaptiveMCMCHandler.cpp.

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

◆ CalculateRobbinsMonroStepLength()

void adaptive_mcmc::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 516 of file AdaptiveMCMCHandler.cpp.

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

◆ CheckMatrixValidityForAdaption()

void adaptive_mcmc::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 450 of file AdaptiveMCMCHandler.cpp.

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

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

Definition at line 144 of file AdaptiveMCMCHandler.cpp.

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

◆ CurrVal()

double adaptive_mcmc::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 508 of file AdaptiveMCMCHandler.cpp.

508  {
509 // ********************************************
512  return (*_fCurrVal)[par_index];
513 }
const std::vector< double > * _fCurrVal
Current values of parameters.

◆ GetAdaptionScale()

double adaptive_mcmc::AdaptiveMCMCHandler::GetAdaptionScale ( )
inline

Get the current adaption scale.

Definition at line 223 of file AdaptiveMCMCHandler.h.

224  {
225  return adaption_scale;
226  }

◆ GetAdaptiveCovariance()

TMatrixDSym* adaptive_mcmc::AdaptiveMCMCHandler::GetAdaptiveCovariance ( ) const
inline

Increase by one number of total steps.

Definition at line 203 of file AdaptiveMCMCHandler.h.

204  {
205  return adaptive_covariance;
206  }

◆ GetNFixed()

int adaptive_mcmc::AdaptiveMCMCHandler::GetNFixed ( ) const
inline

Get the number of fixed parameters.

Definition at line 140 of file AdaptiveMCMCHandler.h.

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

◆ GetNumParams()

int adaptive_mcmc::AdaptiveMCMCHandler::GetNumParams ( ) const
inline

Get the current values of the parameters.

Definition at line 159 of file AdaptiveMCMCHandler.h.

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

◆ GetOutFileName()

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

Get Name of Output File.

Definition at line 217 of file AdaptiveMCMCHandler.h.

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

◆ GetParameterMeans()

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

Get the parameter means used in the adaptive handler.

Definition at line 210 of file AdaptiveMCMCHandler.h.

211  {
212  return par_means;
213  }

◆ GetParIndex()

int adaptive_mcmc::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 97 of file AdaptiveMCMCHandler.h.

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

◆ GetParIndicesFromNames()

std::vector<int> adaptive_mcmc::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 112 of file AdaptiveMCMCHandler.h.

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

◆ GetTotalSteps()

int adaptive_mcmc::AdaptiveMCMCHandler::GetTotalSteps ( ) const
inline

Get Total Number of Steps.

Definition at line 179 of file AdaptiveMCMCHandler.h.

180  {
181  return total_steps;
182  }

◆ GetUseRobbinsMonro()

bool adaptive_mcmc::AdaptiveMCMCHandler::GetUseRobbinsMonro ( ) const
inline

Use Robbins-Monro approach?

Definition at line 229 of file AdaptiveMCMCHandler.h.

230  {
231  return use_robbins_monro;
232  }

◆ IncrementAcceptedSteps()

void adaptive_mcmc::AdaptiveMCMCHandler::IncrementAcceptedSteps ( )
inline

Definition at line 196 of file AdaptiveMCMCHandler.h.

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

◆ IncrementNSteps()

void adaptive_mcmc::AdaptiveMCMCHandler::IncrementNSteps ( )
inline

Increase by one number of total steps.

Definition at line 191 of file AdaptiveMCMCHandler.h.

192  {
193  total_steps++;
194  }

◆ IndivStepScaleAdapt()

bool adaptive_mcmc::AdaptiveMCMCHandler::IndivStepScaleAdapt ( ) const

Tell whether we want reset step scale or not.

Definition at line 401 of file AdaptiveMCMCHandler.cpp.

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

◆ InitFromConfig()

bool adaptive_mcmc::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 33 of file AdaptiveMCMCHandler.cpp.

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

◆ IsFixed()

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

Check if a parameter is fixed.

Definition at line 165 of file AdaptiveMCMCHandler.h.

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

◆ Print()

void adaptive_mcmc::AdaptiveMCMCHandler::Print ( ) const

Print all class members.

Definition at line 435 of file AdaptiveMCMCHandler.cpp.

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

◆ SaveAdaptiveToFile()

void adaptive_mcmc::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 194 of file AdaptiveMCMCHandler.cpp.

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

◆ SetAdaptiveBlocks()

void adaptive_mcmc::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 152 of file AdaptiveMCMCHandler.cpp.

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

◆ SetFixed()

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

Set the fixed parameters.

Definition at line 133 of file AdaptiveMCMCHandler.h.

134  {
135  _fFixedPars = fix;
136  }

◆ SetParams()

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

Set the current values of the parameters.

Definition at line 127 of file AdaptiveMCMCHandler.h.

128  {
129  _fCurrVal = params;
130  }

◆ SetThrowMatrixFromFile()

void adaptive_mcmc::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 280 of file AdaptiveMCMCHandler.cpp.

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

◆ SetTotalSteps()

void adaptive_mcmc::AdaptiveMCMCHandler::SetTotalSteps ( const int  nsteps)
inline

Change Total Number of Steps to new value.

Definition at line 185 of file AdaptiveMCMCHandler.h.

186  {
187  total_steps = nsteps;
188  }

◆ SkipAdaption()

bool adaptive_mcmc::AdaptiveMCMCHandler::SkipAdaption ( ) const

Tell if we are Skipping Adaption.

Definition at line 420 of file AdaptiveMCMCHandler.cpp.

420  {
421 // ********************************************
423  total_steps< start_adaptive_update) return true;
424  else return false;
425 }

◆ UpdateAdaptiveCovariance()

void adaptive_mcmc::AdaptiveMCMCHandler::UpdateAdaptiveCovariance ( )

Method to update adaptive MCMC [13].

Parameters
_fCurrValValue of each parameter necessary for updating throw matrix

Definition at line 332 of file AdaptiveMCMCHandler.cpp.

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

◆ UpdateMatrixAdapt()

bool adaptive_mcmc::AdaptiveMCMCHandler::UpdateMatrixAdapt ( )

Tell whether matrix should be updated.

Definition at line 408 of file AdaptiveMCMCHandler.cpp.

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

◆ UpdateRobbinsMonroScale()

void adaptive_mcmc::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 542 of file AdaptiveMCMCHandler.cpp.

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

Member Data Documentation

◆ _fCurrVal

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

Current values of parameters.

Definition at line 287 of file AdaptiveMCMCHandler.h.

◆ _fFixedPars

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

Vector of fixed parameters.

Definition at line 284 of file AdaptiveMCMCHandler.h.

◆ _param_skip_adapt_flags

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

Parameters to skip during adaption.

Definition at line 314 of file AdaptiveMCMCHandler.h.

◆ _ParamNames

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

Parameter names.

Definition at line 290 of file AdaptiveMCMCHandler.h.

◆ acceptance_rate_batch_size

int adaptive_mcmc::AdaptiveMCMCHandler::acceptance_rate_batch_size
private

Acceptance rate in the current batch.

Definition at line 293 of file AdaptiveMCMCHandler.h.

◆ adapt_block_matrix_indices

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

Indices for block-matrix adaption.

Definition at line 262 of file AdaptiveMCMCHandler.h.

◆ adapt_block_sizes

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

Size of blocks for adaption.

Definition at line 265 of file AdaptiveMCMCHandler.h.

◆ adaption_scale

double adaptive_mcmc::AdaptiveMCMCHandler::adaption_scale
private

Scaling factor.

Definition at line 278 of file AdaptiveMCMCHandler.h.

◆ adaptive_covariance

TMatrixDSym* adaptive_mcmc::AdaptiveMCMCHandler::adaptive_covariance
private

Full adaptive covariance matrix.

Definition at line 272 of file AdaptiveMCMCHandler.h.

◆ adaptive_save_n_iterations

int adaptive_mcmc::AdaptiveMCMCHandler::adaptive_save_n_iterations
private

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

Definition at line 256 of file AdaptiveMCMCHandler.h.

◆ adaptive_update_step

int adaptive_mcmc::AdaptiveMCMCHandler::adaptive_update_step
private

Steps between changing throw matrix.

Definition at line 252 of file AdaptiveMCMCHandler.h.

◆ c_robbins_monro

double adaptive_mcmc::AdaptiveMCMCHandler::c_robbins_monro
private

Constant "step scaling" factor for Robbins-Monro.

Definition at line 302 of file AdaptiveMCMCHandler.h.

◆ end_adaptive_update

int adaptive_mcmc::AdaptiveMCMCHandler::end_adaptive_update
private

Steps between changing throw matrix.

Definition at line 249 of file AdaptiveMCMCHandler.h.

◆ initial_scale

double adaptive_mcmc::AdaptiveMCMCHandler::initial_scale
private

Initial scaling factor.

Definition at line 281 of file AdaptiveMCMCHandler.h.

◆ initial_throw_matrix

const TMatrixDSym* adaptive_mcmc::AdaptiveMCMCHandler::initial_throw_matrix
private

Initial throw matrix.

Definition at line 317 of file AdaptiveMCMCHandler.h.

◆ initial_throw_matrix_saved

bool adaptive_mcmc::AdaptiveMCMCHandler::initial_throw_matrix_saved = false
private

Flag to check if initial throw matrix has been saved.

Definition at line 320 of file AdaptiveMCMCHandler.h.

◆ n_rm_restarts

int adaptive_mcmc::AdaptiveMCMCHandler::n_rm_restarts
private

Number of restarts for Robbins Monro (so far)

Definition at line 305 of file AdaptiveMCMCHandler.h.

◆ output_file_name

std::string adaptive_mcmc::AdaptiveMCMCHandler::output_file_name
private

Name of the file to save the adaptive matrices into.

Definition at line 259 of file AdaptiveMCMCHandler.h.

◆ par_means

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

Mean values for all parameters.

Definition at line 269 of file AdaptiveMCMCHandler.h.

◆ prev_step_accepted

bool adaptive_mcmc::AdaptiveMCMCHandler::prev_step_accepted
private

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

Definition at line 311 of file AdaptiveMCMCHandler.h.

◆ start_adaptive_throw

int adaptive_mcmc::AdaptiveMCMCHandler::start_adaptive_throw
private

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

Definition at line 243 of file AdaptiveMCMCHandler.h.

◆ start_adaptive_update

int adaptive_mcmc::AdaptiveMCMCHandler::start_adaptive_update
private

When do we stop update the adaptive matrix.

Definition at line 246 of file AdaptiveMCMCHandler.h.

◆ target_acceptance

double adaptive_mcmc::AdaptiveMCMCHandler::target_acceptance
private

Target acceptance rate for Robbins Monro.

Definition at line 299 of file AdaptiveMCMCHandler.h.

◆ total_rm_restarts

int adaptive_mcmc::AdaptiveMCMCHandler::total_rm_restarts
private

Total number of restarts ALLOWED for Robbins Monro.

Definition at line 308 of file AdaptiveMCMCHandler.h.

◆ total_steps

int adaptive_mcmc::AdaptiveMCMCHandler::total_steps
private

Total number of MCMC steps.

Definition at line 275 of file AdaptiveMCMCHandler.h.

◆ use_robbins_monro

bool adaptive_mcmc::AdaptiveMCMCHandler::use_robbins_monro
private

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

Definition at line 296 of file AdaptiveMCMCHandler.h.


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