Static Adaptation
Defines
-
__SMC_STATICMODELADAPT_HH
-
namespace smc
The Sequential Monte Carlo namespace.
The classes and functions within this namespace are intended to be used for producing implemenetations of SMC samplers and related simulation techniques.
-
class staticModelAdapt
- #include <staticModelAdapt.h>
A class containing parameters and functions to update these parameters in the context of static Bayesian models.
Public Functions
-
inline ~staticModelAdapt()
-
inline staticModelAdapt()
The class constructor which sets the current and previous temperatures to zero.
-
inline void ChooseTemp(const arma::vec &logweight, const arma::vec &loglike, double desiredCESS, double epsilon = 0.01)
Chooses the next temperature such that a desired conditional ESS is maintained.
Chooses the next temperature such that a desired conditional ESS is maintained.
- Parameters
logweight – An armadillo vector containing the logarithm of the current particle weights.
loglike – An armadillo vector containing the log likelihood of the current particle values.
desiredCESS – The target conditional ESS for the next temperature (generally fixed).
epsilon – The tolerance for the bisection method (maximum difference between desired and actual conditional ESS).
-
inline void calcEmpCov(const arma::mat &theta, const arma::vec &logweight)
Calculates the empirical covariance matrix based on the current weighted particle set.
Calculates the empirical covariance matrix based on the current weighted particle set.
- Parameters
theta – An [Nxd] armadillo matrix of doubles for the current particle values, where N is the number of particles and d is the dimension of the parameter.
logweight – An armadillo vector containing the logarithm of the current particle weights.
-
inline void calcCholCov(const arma::mat &theta, const arma::vec logweight)
Calculates the Cholesky decomposition of the empirical covariance matrix based on the current weighted particle set.
Calculates the Cholesky decomposition of the empirical covariance matrix based on the current weighted particle set.
- Parameters
theta – An [Nxd] armadillo matrix of doubles for the current particle values, where N is the number of particles and d is the dimension of the parameter
logweight – An armadillo vector of the logarithm of the current particle weights
-
inline int calcMcmcRepeats(double acceptProb, double desiredAcceptProb, int initialN, int maxReps)
Calculates the number of MCMC repeats based on the results of the most recent set of MCMC moves.
Calculates the number of MCMC repeats based on the results of the most recent set of MCMC moves.
- Parameters
acceptProb – The proportion of accepted MCMC steps in the most recent iteration.
desiredAcceptProb – The desired probability of a successful move for each particle.
initialN – The initial number of MCMC repeats.
maxReps – The maximum number of MCMC repeats.
-
inline double GetTemp(int t) const
Returns the t-th temperature.
-
inline double GetTempCurr(void) const
Returns the current temperature.
-
inline double GetTempPrevious(void) const
Returns the previous temperature.
-
inline arma::mat GetCholCov(void) const
Returns the Cholesky decomposition of the empirical covariance matrix based on the current weighted particle set.
-
inline arma::mat GetEmpCov(void) const
Returns the empirical covariance matrix based on the current weighted particle set.
Private Functions
-
inline double CESSdiff(const arma::vec &logweight, const arma::vec &loglike, double tempDiff, double desiredCESS)
Computes the difference between the conditional ESS given the specified temperature difference and the desired conditional ESS.
-
inline double bisection(double curr, const arma::vec &logweight, const arma::vec &loglike, double desiredCESS, double epsilon)
Performs the bisection method to find the temperature within (temp_curr,1) which gives the desired conditional ESS.
-
inline ~staticModelAdapt()
-
class staticModelAdapt