normalize_logpost {aghq}R Documentation

Normalize the joint posterior using AGHQ

Description

This function takes in the optimization results from aghq::optimize_theta() and returns a list with the quadrature points, weights, and normalization information. Like aghq::optimize_theta(), this is designed for use only within aghq::aghq, but is exported for debugging and documented in case you want to modify it somehow, or something.

Usage

normalize_logpost(
  optresults,
  k,
  whichfirst = 1,
  basegrid = NULL,
  ndConstruction = "product",
  ...
)

Arguments

optresults

The results of calling aghq::optimize_theta(): see return value of that function.

k

Integer, the number of quadrature points to use. I suggest at least 3. k = 1 corresponds to a Laplace approximation.

whichfirst

Integer between 1 and the dimension of the parameter space, default 1. The user shouldn't have to worry about this: it's used internally to re-order the parameter vector before doing the quadrature, which is useful when calculating marginal posteriors.

basegrid

Optional. Provide an object of class NIGrid from the mvQuad package, representing the base quadrature rule that will be adapted. This is only for users who want more complete control over the quadrature, and is not necessary if you are fine with the default option which basically corresponds to mvQuad::createNIGrid(length(theta),'GHe',k,'product').

ndConstruction

Create a multivariate grid using a product or sparse construction? Passed directly to mvQuad::createNIGrid(), see that function for further details. Note that the use of sparse grids within aghq is currently experimental and not supported by tests. In particular, calculation of marginal posteriors is known to fail currently.

...

Additional arguments to be passed to optresults$ff, see ?optimize_theta.

Value

If k > 1, a list with elements:

See Also

Other quadrature: aghq(), get_hessian(), get_log_normconst(), get_mode(), get_nodesandweights(), get_numquadpoints(), get_opt_results(), get_param_dim(), laplace_approximation(), marginal_laplace_tmb(), marginal_laplace(), nested_quadrature(), optimize_theta(), plot.aghq(), print.aghqsummary(), print.aghq(), print.laplacesummary(), print.laplace(), print.marginallaplacesummary(), summary.aghq(), summary.laplace(), summary.marginallaplace()

Examples

# Same setup as optimize_theta
logfteta <- function(eta,y) {
  sum(y) * eta - (length(y) + 1) * exp(eta) - sum(lgamma(y+1)) + eta
}
set.seed(84343124)
y <- rpois(10,5) # Mode should be sum(y) / (10 + 1)
truemode <- log((sum(y) + 1)/(length(y) + 1))
objfunc <- function(x) logfteta(x,y)
funlist <- list(
  fn = objfunc,
  gr = function(x) numDeriv::grad(objfunc,x),
  he = function(x) numDeriv::hessian(objfunc,x)
)
opt_sparsetrust <- optimize_theta(funlist,1.5)
opt_trust <- optimize_theta(funlist,1.5,control = default_control(method = "trust"))
opt_bfgs <- optimize_theta(funlist,1.5,control = default_control(method = "BFGS"))

# Quadrature with 3, 5, and 7 points using sparse trust region optimization:
norm_sparse_3 <- normalize_logpost(opt_sparsetrust,3,1)
norm_sparse_5 <- normalize_logpost(opt_sparsetrust,5,1)
norm_sparse_7 <- normalize_logpost(opt_sparsetrust,7,1)

# Quadrature with 3, 5, and 7 points using dense trust region optimization:
norm_trust_3 <- normalize_logpost(opt_trust,3,1)
norm_trust_5 <- normalize_logpost(opt_trust,5,1)
norm_trust_7 <- normalize_logpost(opt_trust,7,1)

# Quadrature with 3, 5, and 7 points using BFGS optimization:
norm_bfgs_3 <- normalize_logpost(opt_bfgs,3,1)
norm_bfgs_5 <- normalize_logpost(opt_bfgs,5,1)
norm_bfgs_7 <- normalize_logpost(opt_bfgs,7,1)


[Package aghq version 0.4.1 Index]