Selecting the Box-Cox parameter for general d using Rcpp
Source:R/box_cox_functions_rcpp.R
find_lambda_rcpp.RdFinds a value of the Box-Cox transformation parameter lambda for which
the (positive) random variable with log-density \(\log f\) has a
density closer to that of a Gaussian random variable.
In the following we use theta (\(\theta\)) to denote the argument
of logf on the original scale and phi (\(\phi\)) on the
Box-Cox transformed scale.
Arguments
- logf
A pointer to a compiled C++ function returning the log of the target density \(f\).
- ...
further arguments to be passed to
logfand related functions.- d
A numeric scalar. Dimension of \(f\).
- n_grid
A numeric scalar. Number of ordinates for each variable in
phi. If this is not supplied a default value ofceiling(2501 ^ (1 / d))is used.- ep_bc
A (positive) numeric scalar. Smallest possible value of
phito consider. Used to avoid negative values ofphi.- min_phi, max_phi
Numeric vectors. Smallest and largest values of
phiat which to evaluatelogf, i.e., the range of values ofphiover which to evaluatelogf. Any components inmin_phithat are not positive are set toep_bc.- which_lam
A numeric vector. Contains the indices of the components of
phithat ARE to be Box-Cox transformed.- lambda_range
A numeric vector of length 2. Range of lambda over which to optimise.
- init_lambda
A numeric vector of length 1 or
d. Initial value of lambda used in the search for the best lambda. Ifinit_lambdais a scalar thenrep(init_lambda, d)is used.- phi_to_theta
A pointer to a compiled C++ function returning (the inverse) of the transformation from
thetatophiused to ensure positivity ofphiprior to Box-Cox transformation. The argument isphiand the returned value istheta. Ifphi_to_thetais undefined at the input value then the function should returnNA.- log_j
A pointer to a compiled C++ function returning the log of the Jacobian of the transformation from
thetatophi, i.e., based on derivatives of \(phi\) with respect to \(theta\). Takesthetaas its argument.- user_args
A list of numeric components providing arguments to the user-supplied functions
phi_to_thetaandlog_j.
Value
A list containing the following components
- lambda
A numeric vector. The value of
lambda.- gm
A numeric vector. Box-Cox scaling parameter, estimated by the geometric mean of the values of phi used in the optimisation to find the value of lambda, weighted by the values of f evaluated at phi.
- init_psi
A numeric vector. An initial estimate of the mode of the Box-Cox transformed density
- sd_psi
A numeric vector. Estimates of the marginal standard deviations of the Box-Cox transformed variables.
- phi_to_theta
as detailed above (only if
phi_to_thetais supplied)- log_j
as detailed above (only if
log_jis supplied)- user_args
as detailed above (only if
user_argsis supplied)
Details
The general idea is to evaluate the density \(f\) on a
d-dimensional grid, with n_grid ordinates for each of the
d variables.
We treat each combination of the variables in the grid as a data point
and perform an estimation of the Box-Cox transformation parameter
lambda, in which each data point is weighted by the density
at that point. The vectors min_phi and max_phi define the
limits of the grid and which_lam can be used to specify that only
certain components of phi are to be transformed.
References
Box, G. and Cox, D. R. (1964) An Analysis of Transformations. Journal of the Royal Statistical Society. Series B (Methodological), 26(2), 211-252.
Andrews, D. F. and Gnanadesikan, R. and Warner, J. L. (1971) Transformations of Multivariate Data, Biometrics, 27(4).
Eddelbuettel, D. and Francois, R. (2011). Rcpp: Seamless R and C++ Integration. Journal of Statistical Software, 40(8), 1-18. doi:10.18637/jss.v040.i08
Eddelbuettel, D. (2013). Seamless R and C++ Integration with Rcpp, Springer, New York. ISBN 978-1-4614-6867-7.
See also
ru_rcpp to perform ratio-of-uniforms sampling.
find_lambda_one_d_rcpp to produce (somewhat)
automatically a list for the argument lambda of ru for the
d = 1 case.
Examples
# Log-normal density ===================
# Note: the default value max_phi = 10 is OK here but this will not always
# be the case
ptr_lnorm <- create_xptr("logdlnorm")
mu <- 0
sigma <- 1
lambda <- find_lambda_rcpp(logf = ptr_lnorm, mu = mu, sigma = sigma)
lambda
#> $lambda
#> [1] 0.05408856
#>
#> $gm
#> [1] 0.971952
#>
#> $init_psi
#> [1] -0.05181524
#>
#> $sd_psi
#> Var1
#> 0.8614544
#>
#> $user_args
#> list()
#>
x <- ru_rcpp(logf = ptr_lnorm, mu = mu, sigma = sigma, d = 1, n = 1000,
trans = "BC", lambda = lambda)
# Gamma density ===================
alpha <- 1
# Choose a sensible value of max_phi
max_phi <- qgamma(0.999, shape = alpha)
# [Of course, typically the quantile function won't be available. However,
# In practice the value of lambda chosen is quite insensitive to the choice
# of max_phi, provided that max_phi is not far too large or far too small.]
ptr_gam <- create_xptr("logdgamma")
lambda <- find_lambda_rcpp(logf = ptr_gam, alpha = alpha, max_phi = max_phi)
lambda
#> $lambda
#> [1] 0.2801406
#>
#> $gm
#> [1] 0.5525366
#>
#> $init_psi
#> [1] -0.2060046
#>
#> $sd_psi
#> Var1
#> 0.573372
#>
#> $user_args
#> list()
#>
x <- ru_rcpp(logf = ptr_gam, alpha = alpha, d = 1, n = 1000, trans = "BC",
lambda = lambda)
# \donttest{
# Generalized Pareto posterior distribution ===================
n <- 1000
# 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)
n <- 1000
# Sample on original scale, with no rotation ----------------
ptr_gp <- create_xptr("loggp")
for_ru_rcpp <- c(list(logf = ptr_gp, init = init, d = 2, n = n,
lower = c(0, -Inf)), ss, rotate = FALSE)
x1 <- do.call(ru_rcpp, for_ru_rcpp)
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.1334443 -0.2113298 0.1711807 0
#> b2minus -0.1072281 0.2724910 -0.1805839 0
#> b1plus 0.1633796 0.2878117 -0.1869423 0
#> b2plus 0.1244646 -0.2071443 0.2242449 0
#>
#> estimated probability of acceptance:
#> [1] 0.1566661
#>
#> sample summary
#> V1 V2
#> Min. :0.5833 Min. :-0.8115
#> 1st Qu.:0.8686 1st Qu.:-0.5570
#> Median :0.9529 Median :-0.4940
#> Mean :0.9573 Mean :-0.4917
#> 3rd Qu.:1.0353 3rd Qu.:-0.4331
#> Max. :1.4199 Max. :-0.1532
# Sample on original scale, with rotation ----------------
for_ru_rcpp <- c(list(logf = ptr_gp, init = init, d = 2, n = n,
lower = c(0, -Inf)), ss)
x2 <- do.call(ru_rcpp, for_ru_rcpp)
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.00000000 0
#> b1minus -0.05520268 -0.08353859 0.03131593 0
#> b2minus -0.06863109 0.06634890 -0.11558233 0
#> b1plus 0.11280290 0.24620479 0.11449588 0
#> b2plus 0.07966329 0.12263834 0.14352731 0
#>
#> estimated probability of acceptance:
#> [1] 0.4253509
#>
#> sample summary
#> V1 V2
#> Min. :0.6111 Min. :-0.80195
#> 1st Qu.:0.8697 1st Qu.:-0.55387
#> Median :0.9460 Median :-0.49226
#> Mean :0.9548 Mean :-0.48938
#> 3rd Qu.:1.0286 3rd Qu.:-0.42586
#> Max. :1.4028 Max. :-0.07772
# Sample on Box-Cox transformed scale ----------------
# Find initial estimates for phi = (phi1, phi2),
# where phi1 = sigma
# and phi2 = xi + sigma / max(x),
# and ranges of phi1 and phi2 over over which to evaluate
# the posterior to find a suitable value of lambda.
temp <- do.call(gpd_init, ss)
min_phi <- pmax(0, temp$init_phi - 2 * temp$se_phi)
max_phi <- pmax(0, temp$init_phi + 2 * temp$se_phi)
# Set phi_to_theta() that ensures positivity of phi
# We use phi1 = sigma and phi2 = xi + sigma / max(data)
# Create an external pointer to this C++ function
ptr_phi_to_theta_gp <- create_phi_to_theta_xptr("gp")
# Note: log_j is set to zero by default inside find_lambda_rcpp()
lambda <- find_lambda_rcpp(logf = ptr_gp, ss = ss, d = 2, min_phi = min_phi,
max_phi = max_phi, user_args = list(xm = ss$xm),
phi_to_theta = ptr_phi_to_theta_gp)
lambda
#> $lambda
#> [1] 0.1785430 0.3966356
#>
#> $gm
#> [1] 0.94787641 0.04372419
#>
#> $init_psi
#> [1] -0.04802372 -0.26940988
#>
#> $sd_psi
#> Var1 Var2
#> 0.11025007 0.03112843
#>
#> $phi_to_theta
#> <pointer: 0x0000021bf53f2260>
#>
#> $log_j
#> <pointer: 0x0000021bf53f2180>
#>
#> $user_args
#> $user_args$xm
#> [1] 1.739151
#>
#>
# Sample on Box-Cox transformed, without rotation
x3 <- ru_rcpp(logf = ptr_gp, ss = ss, d = 2, n = n, trans = "BC",
lambda = lambda, rotate = FALSE)
plot(x3, xlab = "sigma", ylab = "xi")
abline(a = 0, b = -1 / ss$xm)
summary(x3)
#> ru bounding box:
#> box vals1 vals2 conv
#> a 1.00000000 0.0000000 0.00000000 0
#> b1minus -0.14471429 -0.2409108 0.04161396 0
#> b2minus -0.04197202 0.0957779 -0.06588483 0
#> b1plus 0.14425105 0.2398479 -0.02841056 0
#> b2plus 0.04418060 -0.1224707 0.07309138 0
#>
#> estimated probability of acceptance:
#> [1] 0.4798464
#>
#> sample summary
#> V1 V2
#> Min. :0.6166 Min. :-0.82146
#> 1st Qu.:0.8816 1st Qu.:-0.55404
#> Median :0.9489 Median :-0.49610
#> Mean :0.9575 Mean :-0.49240
#> 3rd Qu.:1.0303 3rd Qu.:-0.43350
#> Max. :1.4348 Max. :-0.02356
# Sample on Box-Cox transformed, with rotation
x4 <- ru_rcpp(logf = ptr_gp, ss = ss, d = 2, n = n, trans = "BC",
lambda = lambda)
plot(x4, xlab = "sigma", ylab = "xi")
abline(a = 0, b = -1 / ss$xm)
summary(x4)
#> ru bounding box:
#> box vals1 vals2 conv
#> a 1.00000000 0.000000000 0.000000000 0
#> b1minus -0.07103912 -0.114896067 0.009430772 0
#> b2minus -0.07098231 -0.004547423 -0.111423200 0
#> b1plus 0.07811515 0.132373422 0.006738742 0
#> b2plus 0.07471743 -0.004544138 0.123610798 0
#>
#> estimated probability of acceptance:
#> [1] 0.5277045
#>
#> sample summary
#> V1 V2
#> Min. :0.5947 Min. :-0.8036
#> 1st Qu.:0.8594 1st Qu.:-0.5520
#> Median :0.9414 Median :-0.4828
#> Mean :0.9516 Mean :-0.4845
#> 3rd Qu.:1.0308 3rd Qu.:-0.4237
#> Max. :1.4609 Max. :-0.1185
def_par <- graphics::par(no.readonly = TRUE)
par(mfrow = c(2,2), mar = c(4, 4, 1.5, 1))
plot(x1, xlab = "sigma", ylab = "xi", ru_scale = TRUE,
main = "mode relocation")
plot(x2, xlab = "sigma", ylab = "xi", ru_scale = TRUE,
main = "mode relocation and rotation")
plot(x3, xlab = "sigma", ylab = "xi", ru_scale = TRUE,
main = "Box-Cox and mode relocation")
plot(x4, xlab = "sigma", ylab = "xi", ru_scale = TRUE,
main = "Box-Cox, mode relocation and rotation")
graphics::par(def_par)
# }