mcmcMh: Metropolis-Hastings MCMC

View source: R/mcmc.r

mcmcMhR Documentation

Metropolis-Hastings MCMC

Description

Run nIterations of a Metropolis-Hastings MCMC to sample from the target distribution using a gaussian proposal kernel. Two optional optimizations are also implemented: truncated gaussian proposal (to match the support of the target distribution, i.e. boundary of the parameters) and adaptive gaussian proposal (to match the size and the shape of the target distribution).

Usage

mcmcMh(
  target,
  initTheta,
  proposalSd = NULL,
  nIterations,
  covmat = NULL,
  limits = list(lower = NULL, upper = NULL),
  adaptSizeStart = NULL,
  adaptSizeCooling = 0.99,
  adaptShapeStart = NULL,
  adaptShapeStop = NULL,
  printInfoEvery = nIterations/100,
  verbose = FALSE,
  maxScalingSd = 50
)

Arguments

target

R-function that takes a single argument: theta (named numeric vector of parameter values) and returns a list of 2 elements:

  • logDensity the logged value of the target density, evaluated at theta.

  • trace a named numeric vector of values to be printed in the trace data.frame returned by mcmcMh.

initTheta

named vector of initial parameter values to start the chain.

proposalSd

vector of standard deviations. If this is given and covmat is not, a diagonal matrix will be built from this to use as covariance matrix of the multivariate Gaussian proposal distribution. By default, this is set to initTheta/10.

nIterations

number of iterations to run the MCMC chain.

covmat

named numeric covariance matrix of the multivariate Gaussian proposal distribution. Must have named rows and columns with at least all estimated theta. If proposalSd is given, this is ignored.

limits

limits for the - potentially truncated - multi-variate normal proposal distribution of the MCMC. Contains 2 elements:

  • lower named numeric vector. Lower truncation points in each dimension of the Gaussian proposal distribution. By default they are set to -Inf.

  • upper named numeric vector. Upper truncation points in each dimension of the Gaussian proposal distribution. By default they are set to Inf.

adaptSizeStart

number of iterations to run before adapting the size of the proposal covariance matrix (see note below). Set to NULL (default) if size is not to be adapted.

adaptSizeCooling

cooling factor for the scaling factor of the covariance matrix during size adaptation (see note below).

adaptShapeStart

number of accepted jumps before adapting the shape of the proposal covariance matrix (see note below). Set to NULL (default) if shape is not to be adapted

adaptShapeStop

number of iterations to run with adaptations of the shape of the proposal covariance matrix before stopping. Se to NULL (default) if never to stop.

printInfoEvery

frequency of information on the chain: acceptance rate and state of the chain. Default value to nIterations/100. Set to NULL to avoid any info.

verbose

logical. If TRUE, information are printed.

maxScalingSd

numeric. Maximum value for the scaling factor of the covariance matrix. Avoid too high values for the scaling factor, which might happen due to the exponential update scheme. In this case, the covariance matrix becomes too wide and the sampling from the truncated proposal kernel becomes highly inefficient

Value

a list with 3 elements:

  • trace a data.frame. Each row contains a state of the chain (as returned by target, and an extra column for the logDensity).

  • acceptanceRate acceptance rate of the MCMC chain.

  • covmatProposal last covariance matrix used for proposals.

Note

The size of the proposal covariance matrix is adapted using the following formulae:

\Sigma_{n+1}=\sigma_n * \Sigma_n

with \sigma_n=\sigma_{n-1}*exp(\alpha^n*(acc - 0.234)), where \alpha is equal to adaptSizeCooling and acc is the acceptance rate of the chain.

The shape of the proposal covariance matrix is adapted using the following formulae:

\Sigma_{n+1}=2.38^2/d * \Sigma_n

with \Sigma_n the empirical covariance matrix and d is the number of estimated parameters in the model.

References

Roberts GO, Rosenthal JS. Examples of adaptive MCMC. Journal of Computational and Graphical Statistics. Taylor & Francis; 2009;18(2):349-67.


sbfnk/fitR documentation built on July 18, 2023, 3:28 p.m.