Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_base_rng.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_BASE_RNG_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_BASE_RNG_HPP
3
8
9namespace stan {
10namespace math {
11
34template <typename LLFunc, typename LLArgs, typename CovarFun,
35 typename CovarArgs, bool InitTheta, typename RNG,
36 require_t<is_all_arithmetic_scalar<CovarArgs, LLArgs>>* = nullptr>
37inline Eigen::VectorXd laplace_base_rng(
38 LLFunc&& ll_fun, LLArgs&& ll_args, CovarFun&& covariance_function,
39 CovarArgs&& covar_args, const laplace_options<InitTheta>& options, RNG& rng,
40 std::ostream* msgs) {
41 Eigen::MatrixXd covariance_train = stan::math::apply(
42 [msgs, &covariance_function](auto&&... args) {
43 return covariance_function(std::forward<decltype(args)>(args)..., msgs);
44 },
45 std::forward<CovarArgs>(covar_args));
47 ll_fun, std::forward<LLArgs>(ll_args), covariance_train, options, msgs);
48 Eigen::VectorXd mean_train = covariance_train * md_est.theta_grad;
49 if (options.solver == 1 || options.solver == 2) {
50 Eigen::MatrixXd V_dec
51 = md_est.L.template triangularView<Eigen::Lower>().solve(
52 md_est.W_r * covariance_train);
53 Eigen::MatrixXd Sigma = covariance_train - V_dec.transpose() * V_dec;
54 return multi_normal_rng(std::move(mean_train), std::move(Sigma), rng);
55 } else {
56 Eigen::MatrixXd Sigma
57 = covariance_train
58 - covariance_train
59 * (md_est.W_r
60 - md_est.W_r
61 * md_est.LU.solve(covariance_train * md_est.W_r))
62 * covariance_train;
63 return multi_normal_rng(std::move(mean_train), std::move(Sigma), rng);
64 }
65}
66
67} // namespace math
68} // namespace stan
69
70#endif
StdVectorBuilder< true, Eigen::VectorXd, T_loc >::type multi_normal_rng(const T_loc &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
Return a multivariate normal random variate with the given location and covariance using the specifie...
Reference for calculations of marginal and its gradients: Margossian et al (2020),...
auto 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,...
Eigen::VectorXd laplace_base_rng(LLFunc &&ll_fun, LLArgs &&ll_args, CovarFun &&covariance_function, CovarArgs &&covar_args, const laplace_options< InitTheta > &options, RNG &rng, std::ostream *msgs)
In a latent gaussian model,.
constexpr decltype(auto) apply(F &&f, Tuple &&t, PreArgs &&... pre_args)
Definition apply.hpp:51
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...