Uses the generalized ratio-of-uniforms method to simulate from a distribution with log-density \(\log f\) (up to an additive constant). The density \(f\) must be bounded, perhaps after a transformation of variable.
ru(
logf,
...,
n = 1,
d = 1,
init = NULL,
mode = NULL,
trans = c("none", "BC", "user"),
phi_to_theta = NULL,
log_j = NULL,
user_args = list(),
lambda = rep(1L, d),
lambda_tol = 1e-06,
gm = NULL,
rotate = ifelse(d == 1, FALSE, TRUE),
lower = rep(-Inf, d),
upper = rep(Inf, d),
r = 1/2,
ep = 0L,
a_algor = if (d == 1) "nlminb" else "optim",
b_algor = c("nlminb", "optim"),
a_method = c("Nelder-Mead", "BFGS", "CG", "L-BFGS-B", "SANN", "Brent"),
b_method = c("Nelder-Mead", "BFGS", "CG", "L-BFGS-B", "SANN", "Brent"),
a_control = list(),
b_control = list(),
var_names = NULL,
shoof = 0.2
)
A function returning the log of the target density \(f\)
evaluated at its first argument.
This function should return -Inf
when the density is zero.
It is better to use logf =
explicitly, for example,
ru(logf = dnorm, log = TRUE, init = 0.1)
,
to avoid argument matching problems. In contrast,
ru(dnorm, log = TRUE, init = 0.1)
will throw an error because partial matching results in
logf
being matched to log = TRUE
.
Further arguments to be passed to logf
and related
functions.
A non-negative integer scalar. The number of simulated values
required. If n = 0
then no simulation is performed but the
component box
in the returned object gives the ratio-of-uniforms
bounding box that would have been used.
A positive integer scalar. The dimension of \(f\).
A numeric vector of length d
. Initial estimate of the
mode of logf
.
If trans = "BC"
or trans = "user"
this is after
Box-Cox transformation or user-defined transformation, but before
any rotation of axes.
If init
is not supplied then rep(1, d)
is used.
If length(init) = 1
and d > 1
then
init <- rep(init, length.out = d)
is used.
A numeric vector of length d
. The mode of logf
.
If trans = "BC"
or trans = "user"
this is after
Box-Cox transformation or user-defined transformation, but before
any rotation of axes. Only supply mode
if the mode is known: it
will not be checked. If mode
is supplied then init
is
ignored.
A character scalar. trans = "none"
for no
transformation, trans = "BC"
for Box-Cox transformation,
trans = "user"
for a user-defined transformation.
If trans = "user"
then the transformation should be specified
using phi_to_theta
and log_j
and user_args
may be
used to pass arguments to phi_to_theta
and log_j
.
See Details and the Examples.
A function returning (the inverse) of the transformation
from theta
(\(\theta\)) to phi
(\(\phi\)) that may be
used to ensure positivity of \(\phi\) prior to Box-Cox transformation.
The argument is phi
and the returned value is theta
.
If phi_to_theta
is undefined at the input value then the
function should return NA
. See Details.
If lambda$phi_to_theta
(see argument lambda
below) is
supplied then this is used instead of any function supplied via
phi_to_theta
.
A function returning the log of the Jacobian of the
transformation from theta
(\(\theta\)) to phi
(\(\phi\)),
i.e., based on derivatives of \(\phi\) with respect to \(\theta\).
Takes theta
as its argument.
If lambda$log_j
(see argument lambda
below) is
supplied then this is used instead of any function supplied via
log_j
.
A list of numeric components. If trans = "user"
then user_args
is a list providing arguments to the user-supplied
functions phi_to_theta
and log_j
.
Either
A numeric vector. Box-Cox transformation parameters, or
A list with components
A numeric vector. Box-Cox parameters (required).
A numeric vector. Box-Cox scaling parameters (optional).
If supplied this overrides any gm
supplied by the individual
gm
argument described below.
A numeric vector. Initial estimate of mode after Box-Cox transformation (optional).
A numeric vector. Estimates of the marginal standard deviations of the Box-Cox transformed variables (optional).
As above (optional).
As above (optional).
This list may be created using find_lambda_one_d (for d
= 1)
or find_lambda (for any d
).
A numeric scalar. Any values in lambda that are less
than lambda_tol
in magnitude are set to zero.
A numeric vector. Box-Cox scaling parameters (optional). If
lambda$gm
is supplied in input list lambda
then
lambda$gm
is used, not gm
.
A logical scalar. If TRUE (d
> 1 only) use Choleski
rotation. If d = 1 and rotate = TRUE
then rotate will be set to
FALSE with a warning. See Details.
Numeric vectors. Lower/upper bounds on the arguments of
the function after any transformation from theta to phi implied by
the inverse of phi_to_theta
. If rotate = FALSE
these
are used in all of the optimisations used to construct the bounding box.
If rotate = TRUE
then they are use only in the first optimisation
to maximise the target density.`
If trans = "BC"
components of lower
that are negative are
set to zero without warning and the bounds implied after the Box-Cox
transformation are calculated inside ru
.
A numeric scalar. Parameter of generalized ratio-of-uniforms.
A numeric scalar. Controls initial estimates for optimisations
to find the \(b\)-bounding box parameters. The default (ep
= 0)
corresponds to starting at the mode of logf
small positive values
of ep
move the constrained variable slightly away from the mode in
the correct direction. If ep
is negative its absolute value is
used, with no warning given.
Character scalars. Either "nlminb"
or
"optim"
.
Respective optimisation algorithms used to find \(a(r)\) and
(\(b\)i-(r),
\(b\)i+(r)).
Character scalars. Respective methods used by
optim
to find \(a(r)\) and
(\(b\)i-(r),
\(b\)i+(r)).
Only used if optim
is the chosen algorithm. If d
= 1 then
a_method
and b_method
are set to "Brent"
without
warning.
Lists of control arguments to optim
or
nlminb
to find \(a(r)\) and
(\(b\)i-(r),
\(b\)i+(r))
respectively.
A character (or numeric) vector of length d
. Names
to give to the column(s) of the simulated values.
A numeric scalar in [0, 1]. Sometimes a spurious
non-zero convergence indicator is returned from
optim
or nlminb
).
In this event we try to check that a minimum has indeed been found using
different algorithm. shoof
controls the starting value provided
to this algorithm.
If shoof = 0
then we start from the current solution.
If shoof = 1
then we start from the initial estimate provided
to the previous minimisation. Otherwise, shoof
interpolates
between these two extremes, with a value close to zero giving a starting
value that is close to the current solution.
The exception to this is when the initial and current solutions are equal.
Then we start from the current solution multiplied by 1 - shoof
.
An object of class "ru"
is a list containing the following
components:
An n
by d
matrix of simulated values.
A (2 * d
+ 1) by d
+ 2 matrix of
ratio-of-uniforms bounding box information, with row names indicating
the box parameter. The columns contain
values of box parameters.
d
-1)values of variables at which these box parameters are obtained.
d
convergence indicators.
Scaling of f within ru
and relocation of the
mode to the origin means that the first row of box
will always
be c(1, rep(0, d))
.
A numeric scalar. An estimate of the probability of acceptance.
The value of r
.
The value of d
.
A function. logf
supplied by the user, but
with f scaled by the maximum of the target density used in the
ratio-of-uniforms method (i.e. logf_rho
), to avoid numerical
problems in contouring f in plot.ru
when
d = 2
.
A function. The target function actually used in the ratio-of-uniforms algorithm.
An n
by d
matrix of values simulated
from the function used in the ratio-of-uniforms algorithm.
A list of further arguments to logf
.
The estimated mode of the target density f, after any Box-Cox transformation and/or user supplied transformation, but before mode relocation.
An R function that performs the inverse transformation
from the transformed variable \(\rho\), on which the generalised
ratio-of-uniforms method is performed, back to the original variable
\(\theta\). Note: trans_fn
is not
vectorised with respect to \(\rho\).
For information about the generalised ratio-of-uniforms method and
transformations see the
Introducing rust vignette. This can also be accessed using
vignette("rust-a-vignette", package = "rust")
.
If trans = "none"
and rotate = FALSE
then ru
implements the (multivariate) generalized ratio of uniforms method
described in Wakefield, Gelfand and Smith (1991) using a target
density whose mode is relocated to the origin (`mode relocation') in the
hope of increasing efficiency.
If trans = "BC"
then marginal Box-Cox transformations of each of
the d
variables is performed, with parameters supplied in
lambda
. The function phi_to_theta
may be used, if
necessary, to ensure positivity of the variables prior to Box-Cox
transformation.
If trans = "user"
then the function phi_to_theta
enables
the user to specify their own transformation.
In all cases the mode of the target function is relocated to the origin after any user-supplied transformation and/or Box-Cox transformation.
If d
is greater than one and rotate = TRUE
then a rotation
of the variable axes is performed after mode relocation. The
rotation is based on the Choleski decomposition (see chol) of the
estimated Hessian (computed using optimHess
of the negated
log-density after any user-supplied transformation or Box-Cox
transformation. If any of the eigenvalues of the estimated Hessian are
non-positive (which may indicate that the estimated mode of logf
is close to a variable boundary) then rotate
is set to FALSE
with a warning. A warning is also given if this happens when
d
= 1.
The default value of the tuning parameter r
is 1/2, which is
likely to be close to optimal in many cases, particularly if
trans = "BC"
.
Wakefield, J. C., Gelfand, A. E. and Smith, A. F. M. (1991) Efficient generation of random variates via the ratio-of-uniforms method. Statistics and Computing (1991), 1, 129-133. doi:10.1007/BF01889987 .
ru_rcpp
for a version of ru
that uses
the Rcpp package to improve efficiency.
summary.ru
for summaries of the simulated values
and properties of the ratio-of-uniforms algorithm.
plot.ru
for a diagnostic plot.
find_lambda_one_d
to produce (somewhat) automatically
a list for the argument lambda
of ru
for the
d
= 1 case.
find_lambda
to produce (somewhat) automatically
a list for the argument lambda
of ru
for any value of
d
.
optim
for choices of the arguments
a_method
, b_method
, a_control
and b_control
.
nlminb
for choices of the arguments
a_control
and b_control
.
optimHess
for Hessian estimation.
chol
for the Choleski decomposition.
# Normal density ===================
# One-dimensional standard normal ----------------
x <- ru(logf = function(x) -x ^ 2 / 2, d = 1, n = 1000, init = 0.1)
# Two-dimensional standard normal ----------------
x <- ru(logf = function(x) -(x[1]^2 + x[2]^2) / 2, d = 2, n = 1000,
init = c(0, 0))
# Two-dimensional normal with positive association ----------------
rho <- 0.9
covmat <- matrix(c(1, rho, rho, 1), 2, 2)
log_dmvnorm <- function(x, mean = rep(0, d), sigma = diag(d)) {
x <- matrix(x, ncol = length(x))
d <- ncol(x)
- 0.5 * (x - mean) %*% solve(sigma) %*% t(x - mean)
}
# No rotation.
x <- ru(logf = log_dmvnorm, sigma = covmat, d = 2, n = 1000, init = c(0, 0),
rotate = FALSE)
# With rotation.
x <- ru(logf = log_dmvnorm, sigma = covmat, d = 2, n = 1000, init = c(0, 0))
# three-dimensional normal with positive association ----------------
covmat <- matrix(rho, 3, 3) + diag(1 - rho, 3)
# No rotation. Slow !
x <- ru(logf = log_dmvnorm, sigma = covmat, d = 3, n = 1000,
init = c(0, 0, 0), rotate = FALSE)
# With rotation.
x <- ru(logf = log_dmvnorm, sigma = covmat, d = 3, n = 1000,
init = c(0, 0, 0))
# Log-normal density ===================
# Sampling on original scale ----------------
x <- ru(logf = dlnorm, log = TRUE, d = 1, n = 1000, lower = 0, init = 1)
# Box-Cox transform with lambda = 0 ----------------
lambda <- 0
x <- ru(logf = dlnorm, log = TRUE, d = 1, n = 1000, lower = 0, init = 0.1,
trans = "BC", lambda = lambda)
# Equivalently, we could use trans = "user" and supply the (inverse) Box-Cox
# transformation and the log-Jacobian by hand
x <- ru(logf = dlnorm, log = TRUE, d = 1, n = 1000, init = 0.1,
trans = "user", phi_to_theta = function(x) exp(x),
log_j = function(x) -log(x))
# Gamma(alpha, 1) density ===================
# Note: the gamma density in unbounded when its shape parameter is < 1.
# Therefore, we can only use trans="none" if the shape parameter is >= 1.
# Sampling on original scale ----------------
alpha <- 10
x <- ru(logf = dgamma, shape = alpha, log = TRUE, d = 1, n = 1000,
lower = 0, init = alpha)
alpha <- 1
x <- ru(logf = dgamma, shape = alpha, log = TRUE, d = 1, n = 1000,
lower = 0, init = alpha)
#> Warning: The Hessian of the target log-density at its mode is not positive
#> definite. This may not be a problem, but it may be that a mode
#> at/near a parameter boundary has been found and/or that the target
#> function is unbounded.
#> It might be worth using the option trans = ``BC''.
# Box-Cox transform with lambda = 1/3 works well for shape >= 1. -----------
alpha <- 1
x <- ru(logf = dgamma, shape = alpha, log = TRUE, d = 1, n = 1000,
trans = "BC", lambda = 1/3, init = alpha)
summary(x)
#> ru bounding box:
#> box vals1 conv
#> a 1.000000 0.000000 0
#> b1minus -1.051825 -1.609437 0
#> b1plus 1.096590 1.774103 0
#>
#> estimated probability of acceptance:
#> [1] 0.785546
#>
#> sample summary
#> V1
#> Min. :0.000494
#> 1st Qu.:0.289279
#> Median :0.695389
#> Mean :1.000469
#> 3rd Qu.:1.387574
#> Max. :6.480577
# Equivalently, we could use trans = "user" and supply the (inverse) Box-Cox
# transformation and the log-Jacobian by hand
# Note: when phi_to_theta is undefined at x this function returns NA
phi_to_theta <- function(x, lambda) {
ifelse(x * lambda + 1 > 0, (x * lambda + 1) ^ (1 / lambda), NA)
}
log_j <- function(x, lambda) (lambda - 1) * log(x)
lambda <- 1/3
x <- ru(logf = dgamma, shape = alpha, log = TRUE, d = 1, n = 1000,
trans = "user", phi_to_theta = phi_to_theta, log_j = log_j,
user_args = list(lambda = lambda), init = alpha)
summary(x)
#> ru bounding box:
#> box vals1 conv
#> a 1.000000 0.000000 0
#> b1minus -1.051825 -1.609437 0
#> b1plus 1.096590 1.774103 0
#>
#> estimated probability of acceptance:
#> [1] 0.7980846
#>
#> sample summary
#> V1
#> Min. :0.001001
#> 1st Qu.:0.302703
#> Median :0.739617
#> Mean :1.049657
#> 3rd Qu.:1.425455
#> Max. :9.801844
# \donttest{
# Generalized Pareto posterior distribution ===================
# Sample data from a GP(sigma, xi) distribution
gpd_data <- rgpd(m = 100, xi = -0.5, sigma = 1)
# Calculate summary statistics for use in the log-likelihood
ss <- gpd_sum_stats(gpd_data)
# Calculate an initial estimate
init <- c(mean(gpd_data), 0)
# Mode relocation only ----------------
n <- 1000
x1 <- ru(logf = gpd_logpost, ss = ss, d = 2, n = n, init = init,
lower = c(0, -Inf), rotate = FALSE)
plot(x1, xlab = "sigma", ylab = "xi")
# Parameter constraint line xi > -sigma/max(data)
# [This may not appear if the sample is far from the constraint.]
abline(a = 0, b = -1 / ss$xm)
summary(x1)
#> ru bounding box:
#> box vals1 vals2 conv
#> a 1.0000000 0.0000000 0.0000000 0
#> b1minus -0.1331405 -0.2116871 0.1524397 0
#> b2minus -0.1010805 0.2842293 -0.1749433 0
#> b1plus 0.1650133 0.2927237 -0.1783390 0
#> b2plus 0.1090747 -0.2109005 0.1991069 0
#>
#> estimated probability of acceptance:
#> [1] 0.1316656
#>
#> sample summary
#> V1 V2
#> Min. :0.6692 Min. :-0.7909
#> 1st Qu.:0.9521 1st Qu.:-0.6055
#> Median :1.0293 Median :-0.5478
#> Mean :1.0378 Mean :-0.5484
#> 3rd Qu.:1.1122 3rd Qu.:-0.4920
#> Max. :1.4167 Max. :-0.2628
# Rotation of axes plus mode relocation ----------------
x2 <- ru(logf = gpd_logpost, ss = ss, d = 2, n = n, init = init,
lower = c(0, -Inf))
plot(x2, xlab = "sigma", ylab = "xi")
abline(a = 0, b = -1 / ss$xm)
summary(x2)
#> ru bounding box:
#> box vals1 vals2 conv
#> a 1.00000000 0.00000000 0.0000000 0
#> b1minus -0.04203948 -0.06214065 0.0359629 0
#> b2minus -0.06194203 0.04729381 -0.1072050 0
#> b1plus 0.11573546 0.26655126 0.1065524 0
#> b2plus 0.06684085 0.12990030 0.1220125 0
#>
#> estimated probability of acceptance:
#> [1] 0.3913894
#>
#> sample summary
#> V1 V2
#> Min. :0.6840 Min. :-0.7967
#> 1st Qu.:0.9523 1st Qu.:-0.6040
#> Median :1.0285 Median :-0.5497
#> Mean :1.0346 Mean :-0.5470
#> 3rd Qu.:1.1144 3rd Qu.:-0.4928
#> Max. :1.4894 Max. :-0.2119
# Cauchy ========================
# The bounding box cannot be constructed if r < 1. For r = 1 the
# bounding box parameters b1-(r) and b1+(r) are attained in the limits
# as x decreases/increases to infinity respectively. This is fine in
# theory but using r > 1 avoids this problem and the largest probability
# of acceptance is obtained for r approximately equal to 1.26.
res <- ru(logf = dcauchy, log = TRUE, init = 0, r = 1.26, n = 1000)
# Half-Cauchy ===================
log_halfcauchy <- function(x) {
return(ifelse(x < 0, -Inf, dcauchy(x, log = TRUE)))
}
# Like the Cauchy case the bounding box cannot be constructed if r < 1.
# We could use r > 1 but the mode is on the edge of the support of the
# density so as an alternative we use a log transformation.
x <- ru(logf = log_halfcauchy, init = 0, trans = "BC", lambda = 0, n = 1000)
x$pa
#> [1] 0.7462687
plot(x, ru_scale = TRUE)
# Example 4 from Wakefield et al. (1991) ===================
# Bivariate normal x bivariate student-t
log_norm_t <- function(x, mean = rep(0, d), sigma1 = diag(d), sigma2 = diag(d)) {
x <- matrix(x, ncol = length(x))
log_h1 <- -0.5 * (x - mean) %*% solve(sigma1) %*% t(x - mean)
log_h2 <- -2 * log(1 + 0.5 * x %*% solve(sigma2) %*% t(x))
return(log_h1 + log_h2)
}
rho <- 0.9
covmat <- matrix(c(1, rho, rho, 1), 2, 2)
y <- c(0, 0)
# Case in the top right corner of Table 3
x <- ru(logf = log_norm_t, mean = y, sigma1 = covmat, sigma2 = covmat,
d = 2, n = 10000, init = y, rotate = FALSE)
x$pa
#> [1] 0.2281282
# Rotation increases the probability of acceptance
x <- ru(logf = log_norm_t, mean = y, sigma1 = covmat, sigma2 = covmat,
d = 2, n = 10000, init = y, rotate = TRUE)
x$pa
#> [1] 0.5216212
# Normal x log-normal: different Box-Cox parameters ==================
norm_lognorm <- function(x, ...) {
dnorm(x[1], ...) + dlnorm(x[2], ...)
}
x <- ru(logf = norm_lognorm, log = TRUE, n = 1000, d = 2, init = c(-1, 0),
trans = "BC", lambda = c(1, 0))
plot(x)
plot(x, ru_scale = TRUE)
# }