Automatic Differentiation
 
Loading...
Searching...
No Matches

◆ laplace_marginal_density_est()

template<typename LLFun , typename LLTupleArgs , typename CovarMat , bool InitTheta, require_t< is_all_arithmetic_scalar< CovarMat, LLTupleArgs > > * = nullptr>
auto stan::math::internal::laplace_marginal_density_est ( LLFun &&  ll_fun,
LLTupleArgs &&  ll_args,
CovarMat &&  covariance,
const laplace_options< InitTheta > &  options,
std::ostream *  msgs 
)
inline

For a latent Gaussian model with hyperparameters phi and latent variables theta, and observations y, this function computes an approximation of the log marginal density, p(y | phi).

This is done by marginalizing out theta, using a Laplace approximation. The latter is obtained by finding the mode, via Newton's method, and computing the Hessian of the likelihood.

The convergence criterion for the Newton/Wolfe loop is a small change in the optimization objective (stored in the Wolfe step state). The user controls the tolerance (i.e. threshold under which the change is deemed small enough) and maximum number of steps.

A description of this algorithm can be found in:

  • (2023) Margossian, "General Adjoint-Differentiated Laplace approximation", https://arxiv.org/pdf/2306.14976. Additional references include:
  • (2020) Margossian et al, "HMC using an adjoint-differentiated Laplace...", NeurIPS, https://arxiv.org/abs/2004.12550.
  • (2006) Rasmussen and Williams, "Gaussian Processes for Machine Learning", second edition, MIT Press, algorithm 3.1.

This function returns intermediate quantities needed by laplace_marginal_density (including its reverse-mode implementation) to compute gradients and/or generate quantities. Which fields are populated depends on the solver that converged (see solver_used).

Template Parameters
LLFunType with a valid operator(ThetaVec, InnerLLTupleArgs) where InnerLLTupleArgs are the elements of LLTupleArgs
LLTupleArgsA tuple whose elements follow the types required for LLFun
CovarFunA functor with an operator()(CovarArgsElements..., {TrainTupleElements...| PredTupleElements...}) method. The operator() method should accept as arguments the inner elements of CovarArgs. The return type of the operator() method should be a type inheriting from Eigen::EigenBase with dynamic sized rows and columns.
CovarArgsA tuple of types to passed as the first arguments of CovarFun::operator()
OpsTupleA tuple of laplace_options types
Parameters
[in]ll_funA log likelihood functor
[in]ll_argsTuple containing parameters for LLFun
[in]covarianceThe covariance matrix for the latent Gaussian
[in]covariance_functiona function which returns the prior covariance.
[in]covar_argsarguments for the covariance function.
[in]optionsA set of options for tuning the solver
[in,out]msgsstream for messages from likelihood and covariance
Returns
A laplace_density_estimates containing:
  1. lmd: the Laplace approximation of the log marginal density, p(y | phi)
  2. theta: the mode of the latent variables
  3. a: the mode in the a parameterization (theta = covariance * a)
  4. theta_grad: gradient of the log-likelihood w.r.t. theta at the mode
  5. W_r: solver-dependent Hessian quantity (see struct docs)
  6. L: solver-dependent factorization of B (Cholesky for solvers 1/2)
  7. LU: LU factorization of B (solver 3 only)
  8. K_root: Cholesky factor of the covariance (solver 2 only)
  9. solver_used: the solver policy that produced the result

Definition at line 1147 of file laplace_marginal_density_estimator.hpp.