Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_marginal.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_PROB_LAPLACE_MARGINAL_HPP
2#define STAN_MATH_MIX_PROB_LAPLACE_MARGINAL_HPP
3
6
7namespace stan {
8namespace math {
30template <bool propto = false, typename LFun, typename LArgs, typename CovarFun,
31 typename CovarArgs, typename OpsTuple>
32inline auto laplace_marginal_tol(LFun&& L_f, LArgs&& l_args,
33 int hessian_block_size,
34 CovarFun&& covariance_function,
35 CovarArgs&& covar_args, OpsTuple&& ops,
36 std::ostream* msgs) {
37 auto options
38 = internal::tuple_to_laplace_options(std::forward<OpsTuple>(ops));
39 options.hessian_block_size = hessian_block_size;
41 std::forward<LFun>(L_f), std::forward<LArgs>(l_args),
42 std::forward<CovarFun>(covariance_function),
43 std::forward<CovarArgs>(covar_args), std::move(options), msgs);
44}
45
65template <bool propto = false, typename LFun, typename LArgs, typename CovarFun,
66 typename CovarArgs>
67inline auto laplace_marginal(LFun&& L_f, LArgs&& l_args, int hessian_block_size,
68 CovarFun&& covariance_function,
69 CovarArgs&& covar_args, std::ostream* msgs) {
70 auto options = laplace_options_default{hessian_block_size};
72 std::forward<LFun>(L_f), std::forward<LArgs>(l_args),
73 std::forward<CovarFun>(covariance_function),
74 std::forward<CovarArgs>(covar_args), options, msgs);
75}
76
77} // namespace math
78} // namespace stan
79
80#endif
Reference for calculations of marginal and its gradients: Margossian et al (2020),...
constexpr auto tuple_to_laplace_options(Options &&ops)
auto laplace_marginal(LFun &&L_f, LArgs &&l_args, int hessian_block_size, CovarFun &&covariance_function, CovarArgs &&covar_args, std::ostream *msgs)
Wrapper function around the laplace_marginal function.
auto laplace_marginal_density(LLFun &&ll_fun, LLTupleArgs &&ll_args, CovarFun &&covariance_function, CovarArgs &&covar_args, const laplace_options< InitTheta > &options, std::ostream *msgs)
For a latent Gaussian model with global parameters phi, latent variables theta, and observations y,...
auto laplace_marginal_tol(LFun &&L_f, LArgs &&l_args, int hessian_block_size, CovarFun &&covariance_function, CovarArgs &&covar_args, OpsTuple &&ops, std::ostream *msgs)
Wrapper function around the laplace_marginal_density function.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...