Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_marginal_density_estimator.hpp File Reference

Detailed Description

Reference for calculations of marginal and its gradients: Margossian et al (2020), https://arxiv.org/abs/2004.12550 and Margossian (2023), https://arxiv.org/pdf/2306.14976.

Definition in file laplace_marginal_density_estimator.hpp.

Go to the source code of this file.

Classes

struct  stan::math::laplace_options_base
 Options for the Laplace approximation. More...
 
struct  stan::math::laplace_options< false >
 
struct  stan::math::laplace_options< true >
 
struct  stan::math::internal::laplace_density_estimates< ThetaVec, WR, L_t, A_vec, ThetaGrad, LU_t, KRoot >
 
struct  stan::math::internal::NewtonState
 Holds the state for the Newton-Raphson optimization loop. More...
 
struct  stan::math::internal::CholeskyWSolverDiag
 Solver Policy 1 (Diagonal): Cholesky decomposition using W. More...
 
struct  stan::math::internal::CholeskyWSolverBlock
 Solver Policy 1 (Block): Cholesky decomposition using block W. More...
 
struct  stan::math::internal::CholeskyKSolver
 Solver Policy 2: Cholesky decomposition of K (Covariance). More...
 
struct  stan::math::internal::LUSolver
 Solver Policy 3: LU Decomposition. More...
 

Namespaces

namespace  stan
 The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation from C or the boost::math::lgamma implementation.
 
namespace  stan::math
 Matrices and templated mathematical functions.
 
namespace  stan::math::internal
 A comparator that works for any container type that has the brackets operator.
 

Typedefs

using stan::math::laplace_options_default = laplace_options< false >
 
using stan::math::laplace_options_user_supplied = laplace_options< true >
 

Functions

auto stan::math::generate_laplace_options (int theta_0_size)
 User function for generating laplace options tuple.
 
template<typename Options >
constexpr auto stan::math::internal::tuple_to_laplace_options (Options &&ops)
 
template<typename WRootMat >
void stan::math::internal::block_matrix_sqrt (WRootMat &W_root, const Eigen::SparseMatrix< double > &W, const Eigen::Index block_size)
 Returns the principal square root of a block diagonal matrix.
 
template<bool InitTheta, typename CovarMat >
void stan::math::internal::validate_laplace_options (const char *frame_name, const laplace_options< InitTheta > &options, const CovarMat &covariance)
 Validates the options for the Laplace approximation.
 
template<typename LLT , typename B_t >
void stan::math::internal::llt_with_jitter (LLT &llt_B, B_t &B, double min_jitter=1e-10, double max_jitter=1e-5)
 Factorize B with jittering fallback.
 
template<typename SolverPolicy , typename NewtonStateT , typename OptionsT , typename LLFunT , typename LLTupleArgsT , typename CovarMatT , typename UpdateFun >
auto stan::math::internal::run_newton_loop (SolverPolicy &solver, NewtonStateT &state, const OptionsT &options, Eigen::Index &step_iter, const LLFunT &ll_fun, const LLTupleArgsT &ll_args, const CovarMatT &covariance, UpdateFun &&update_fun, std::ostream *msgs)
 Run a Newton loop with a solver policy, updating the shared state.
 
void stan::math::internal::log_solver_fallback (const bool allow_fallthrough, std::ostream *msgs, std::string_view context, Eigen::Index iter, std::string_view failed_solver, std::string_view next_solver, const std::exception &e)
 Log a solver fallback event to the provided stream.
 
template<bool InitTheta, typename Opts >
decltype(auto) stan::math::internal::theta_init_impl (Eigen::Index theta_size, Opts &&options)
 
template<typename ObjFun , typename ThetaGradFun , typename Covariance , typename Options >
auto stan::math::internal::create_update_fun (ObjFun &&obj_fun, ThetaGradFun &&theta_grad_f, Covariance &&covariance, Options &&options)
 Create the update function for the line search, capturing necessary references.
 
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)
 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).
 

Variables

static thread_local std::once_flag stan::math::internal::fallback_warning