Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_marginal_density.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_HPP
13#include <unsupported/Eigen/MatrixFunctions>
14#include <cmath>
15
23namespace stan {
24namespace math {
25
53template <
54 typename LLFun, typename LLTupleArgs, typename CovarFun, typename CovarArgs,
55 bool InitTheta,
56 require_t<is_all_arithmetic_scalar<CovarArgs, LLTupleArgs>>* = nullptr>
57inline auto laplace_marginal_density(LLFun&& ll_fun, LLTupleArgs&& ll_args,
58 CovarFun&& covariance_function,
59 CovarArgs&& covar_args,
60 const laplace_options<InitTheta>& options,
61 std::ostream* msgs) {
62 Eigen::MatrixXd covariance = stan::math::apply(
63 [msgs, &covariance_function](auto&&... args) {
64 return covariance_function(std::forward<decltype(args)>(args)..., msgs);
65 },
66 std::forward<CovarArgs>(covar_args));
68 std::forward<LLFun>(ll_fun), std::forward<LLTupleArgs>(ll_args),
69 std::move(covariance), options, msgs)
70 .lmd;
71}
72
73namespace internal {
74
75template <bool ZeroInput = false, typename Output, typename Input,
76 require_tuple_t<Output>* = nullptr, require_tuple_t<Input>* = nullptr,
77 require_t<std::bool_constant<
78 std::tuple_size_v<std::decay_t<Output>> == 0>>* = nullptr,
79 require_t<std::bool_constant<
80 std::tuple_size_v<std::decay_t<Input>> == 0>>* = nullptr>
81inline constexpr void copy_compute_s2(Output&& output, Input&& input) {}
82
91template <bool ZeroInput = false, typename Output, typename Input,
94inline constexpr void copy_compute_s2(Output&& output, Input&& input) {
95 if constexpr (is_tuple_v<Output> && is_tuple_v<Input>) {
96 static_assert(
97 std::tuple_size<std::decay_t<Output>>::value
98 == std::tuple_size<std::decay_t<Input>>::value,
99 "INTERNAL ERROR:(laplace_marginal_lpdf) copy_compute_s2 called on "
100 "tuples of different sizes. This is an internal error, please report "
101 "it: "
102 "https://github.com/stan-dev/math/issues");
103 }
104 return iter_tuple_nested(
105 [](auto&& output_i, auto&& input_i) {
106 using output_i_t = std::decay_t<decltype(output_i)>;
107 if constexpr (is_std_vector_v<output_i_t>) {
108 using dbl_map_t = Eigen::Map<Eigen::Matrix<double, -1, 1>>;
109 using var_map_t = Eigen::Map<Eigen::Matrix<var, -1, 1>>;
110 var_map_t input_map(input_i.data(), input_i.size());
111 dbl_map_t(output_i.data(), output_i.size()).array()
112 += 0.5 * input_map.adj().array();
113 if constexpr (ZeroInput) {
114 input_map.adj().setZero();
115 }
116 } else if constexpr (is_eigen_v<output_i_t>) {
117 output_i.array() += 0.5 * input_i.adj().array();
118 if constexpr (ZeroInput) {
119 input_i.adj().setZero();
120 }
121 } else if constexpr (is_stan_scalar_v<output_i_t>) {
122 output_i += (0.5 * input_i.adj());
123 if constexpr (ZeroInput) {
124 input_i.adj() = 0;
125 }
126 } else {
127 static_assert(
128 sizeof(std::decay_t<output_i_t>*) == 0,
129 "INTERNAL ERROR:(laplace_marginal_lpdf) copy_compute_s2 was "
130 "not able to deduce the actions needed for the given type. "
131 "This is an internal error, please report it: "
132 "https://github.com/stan-dev/math/issues");
133 }
134 },
135 std::forward<Output>(output), std::forward<Input>(input));
136}
137
138} // namespace internal
139
167template <typename LLFun, typename LLTupleArgs, typename CovarFun,
168 typename CovarArgs, bool InitTheta,
170inline auto laplace_marginal_density(LLFun&& ll_fun, LLTupleArgs&& ll_args,
171 CovarFun&& covariance_function,
172 CovarArgs&& covar_args,
173 const laplace_options<InitTheta>& options,
174 std::ostream* msgs) {
175 auto covar_args_refs = to_ref(std::forward<CovarArgs>(covar_args));
176 auto ll_args_refs = to_ref(std::forward<LLTupleArgs>(ll_args));
177 // Solver 1, 2, 3
178 constexpr bool ll_args_contain_var = is_any_var_scalar<LLTupleArgs>::value;
179 auto partial_parm = internal::make_zeroed_arena(ll_args_refs);
180 auto covar_args_adj = internal::make_zeroed_arena(covar_args_refs);
181 double lmd = 0.0;
182 {
183 nested_rev_autodiff nested;
184 // Make one hard copy here
185 auto ll_args_copy = internal::deep_copy_vargs<var>(ll_args_refs);
186 auto covar_args_copy = internal::deep_copy_vargs<var>(covar_args_refs);
187 auto covariance = stan::math::apply(
188 [&covariance_function, &msgs](auto&&... args) {
189 if constexpr (is_any_var_scalar_v<decltype(args)...>) {
190 return to_var_value(covariance_function(args..., msgs));
191 } else {
192 return covariance_function(args..., msgs);
193 }
194 },
195 covar_args_copy);
196 decltype(auto) covariance_val = value_of(covariance);
197 decltype(auto) ll_args_vals = value_of(ll_args_copy);
199 ll_fun, ll_args_vals, covariance_val, options, msgs);
200 auto ll_args_filter = internal::filter_var_scalar_types(ll_args_copy);
201 // tuple of references to var types
202 // Solver 1, 2
203 const bool solver_1_or_2
204 = md_est.solver_used == 1 || md_est.solver_used == 2;
205 arena_t<Eigen::MatrixXd> R(md_est.theta.size() * solver_1_or_2,
206 md_est.theta.size() * solver_1_or_2);
207 // Solver 3
208 arena_t<Eigen::MatrixXd> LU_solve_covariance(
209 covariance.rows() * (md_est.solver_used == 3),
210 covariance.cols() * (md_est.solver_used == 3));
211 // Solver 1, 2, 3
212 arena_t<Eigen::VectorXd> s2(md_est.theta.size());
214 if (md_est.solver_used == 1) {
215 if (options.hessian_block_size == 1) {
216 arena_t<Eigen::MatrixXd> tmp = md_est.W_r.toDense();
217 md_est.L.template triangularView<Eigen::Lower>().solveInPlace(tmp);
218 R.noalias() = tmp.transpose() * tmp;
219 arena_t<Eigen::MatrixXd> C
220 = md_est.L.template triangularView<Eigen::Lower>().solve(
221 md_est.W_r * covariance_val);
222 if constexpr (ll_args_contain_var) {
223 arena_t<Eigen::MatrixXd> A = covariance_val - C.transpose() * C;
224 auto s2_tmp = laplace_likelihood::compute_s2(
225 ll_fun, md_est.theta, A, options.hessian_block_size, ll_args_copy,
226 msgs);
227 s2.deep_copy(s2_tmp);
228 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
229 } else {
230 s2.deep_copy(
231 (0.5
232 * (covariance_val.diagonal() - (C.transpose() * C).diagonal())
233 .cwiseProduct(laplace_likelihood::third_diff(
234 ll_fun, md_est.theta, ll_args_vals, msgs))));
235 }
236
237 } else {
238 arena_t<Eigen::MatrixXd> tmp = md_est.W_r.toDense();
239 md_est.L.template triangularView<Eigen::Lower>().solveInPlace(tmp);
240 R.noalias() = tmp.transpose() * tmp;
241 arena_t<Eigen::MatrixXd> C
242 = md_est.L.template triangularView<Eigen::Lower>().solve(
243 md_est.W_r * covariance_val);
244 arena_t<Eigen::MatrixXd> A = covariance_val - C.transpose() * C;
245 auto s2_tmp = laplace_likelihood::compute_s2(ll_fun, md_est.theta, A,
246 options.hessian_block_size,
247 ll_args_copy, msgs);
248 s2.deep_copy(s2_tmp);
249 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
250 }
251 } else if (md_est.solver_used == 2) {
252 R = md_est.W_r
253 - md_est.W_r * md_est.K_root
254 * md_est.L.transpose()
255 .template triangularView<Eigen::Upper>()
256 .solve(
257 md_est.L.template triangularView<Eigen::Lower>()
258 .solve(md_est.K_root.transpose() * md_est.W_r));
259
260 arena_t<Eigen::MatrixXd> C
261 = md_est.L.template triangularView<Eigen::Lower>().solve(
262 md_est.K_root.transpose());
263 auto s2_tmp = laplace_likelihood::compute_s2(
264 ll_fun, md_est.theta, (C.transpose() * C).eval(),
265 options.hessian_block_size, ll_args_copy, msgs);
266 s2.deep_copy(s2_tmp);
267 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
268 } else { // options.solver with LU decomposition
269 LU_solve_covariance = md_est.LU.solve(covariance_val);
270 auto I_minus_BinvKW
271 = Eigen::MatrixXd::Identity(md_est.W_r.rows(), md_est.W_r.cols())
272 - LU_solve_covariance * md_est.W_r;
273 R = md_est.W_r * I_minus_BinvKW; // == W - W B^{-1} K W
274 arena_t<Eigen::MatrixXd> A
275 = covariance_val - covariance_val * md_est.W_r * LU_solve_covariance;
276 auto s2_tmp = laplace_likelihood::compute_s2(ll_fun, md_est.theta, A,
277 options.hessian_block_size,
278 ll_args_copy, msgs);
279 s2.deep_copy(s2_tmp);
280 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
281 }
282 if constexpr (is_any_var_scalar_v<scalar_type_t<CovarArgs>>) {
283 arena_t<Eigen::MatrixXd> K_adj_arena
284 = 0.5 * md_est.a * md_est.a.transpose() - 0.5 * R
285 + s2 * md_est.theta_grad.transpose()
286 - (R * (covariance.val() * s2)) * md_est.theta_grad.transpose();
288 0.0, [covariance, K_adj_arena](auto&& vi) mutable {
289 covariance.adj().array() += vi.adj() * K_adj_arena.array();
290 });
291 grad(Z.vi_);
292 auto covar_args_filter
293 = internal::filter_var_scalar_types(covar_args_copy);
294 internal::collect_adjoints(covar_args_adj, covar_args_filter);
295 }
296 if constexpr (ll_args_contain_var) {
297 laplace_likelihood::ll_arg_grad(ll_fun, md_est.theta, ll_args_copy, msgs);
298 internal::collect_adjoints<ZeroOut>(partial_parm, ll_args_filter);
299 arena_t<Eigen::VectorXd> v;
300 if (md_est.solver_used == 1 || md_est.solver_used == 2) {
301 v = covariance_val * s2 - covariance_val * R * covariance_val * s2;
302 } else {
303 v = LU_solve_covariance * s2;
304 }
305 laplace_likelihood::diff_eta_implicit(ll_fun, v, md_est.theta,
306 ll_args_copy, msgs);
307 internal::collect_adjoints<ZeroOut>(partial_parm, ll_args_filter);
308 }
309 lmd = md_est.lmd;
310 }
311 var ret(lmd);
312 if constexpr (is_any_var_scalar_v<CovarArgs>) {
313 auto covar_args_filter = internal::filter_var_scalar_types(covar_args_refs);
314 internal::reverse_pass_collect_adjoints(ret, covar_args_filter,
315 std::move(covar_args_adj));
316 }
317 if constexpr (ll_args_contain_var) {
318 auto ll_args_filter = internal::filter_var_scalar_types(ll_args_refs);
319 internal::reverse_pass_collect_adjoints(ret, ll_args_filter,
320 std::move(partial_parm));
321 }
322 return ret;
323}
324
325} // namespace math
326} // namespace stan
327
328#endif
Reference for calculations of marginal and its gradients: Margossian et al (2020),...
constexpr decltype(auto) filter_var_scalar_types(T &&t)
Filter a tuple and return a tuple with references to the types with a var scalar type.
void reverse_pass_collect_adjoints(var ret, Output &&output, Input &&input)
Collects adjoints from a tuple or std::vector of tuples.
constexpr auto make_zeroed_arena(Input &&input)
Creates an arena type that is the same type as the input and initialized with zeros.
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,...
void collect_adjoints(Output &output, Input &&input)
Collect the adjoints from the input and add them to the output.
constexpr void copy_compute_s2(Output &&output, Input &&input)
auto compute_s2(F &&f, Theta &&theta, AMat &&A, int hessian_block_size, TupleArgs &&ll_args, Stream *msgs)
A wrapper that accepts a tuple as arguments.
auto diff_eta_implicit(F &&f, V_t &&v, Theta &&theta, TupleArgs &&ll_args, Stream *msgs)
A wrapper that accepts a tuple as arguments.
auto ll_arg_grad(F &&f, Theta &&theta, TupleArgs &&ll_tup, Stream *msgs=nullptr)
Eigen::VectorXd third_diff(F &&f, Theta &&theta, TupleArgs &&ll_args, Stream *msgs)
A wrapper that accepts a tuple as arguments.
void iter_tuple_nested(F &&f, Types &&... args)
Iterate and nest into a tuple or std::vector to apply f to each matrix or scalar type.
var_value< plain_type_t< T > > make_callback_var(T &&value, F &&functor)
Creates a new var initialized with a callback_vari with a given value and reverse-pass callback funct...
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition value_of.hpp:18
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,...
var_value< Eigen::Matrix< double, T::RowsAtCompileTime, T::ColsAtCompileTime > > to_var_value(const T &a)
Converts an Eigen matrix (or vector or row_vector) or expression of vars into var_value.
var_value< double > var
Definition var.hpp:1187
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:18
static void grad()
Compute the gradient for all variables starting from the end of the AD tape.
Definition grad.hpp:26
constexpr decltype(auto) apply(F &&f, Tuple &&t, PreArgs &&... pre_args)
Definition apply.hpp:51
constexpr bool is_any_var_scalar_v
Definition is_var.hpp:50
std::enable_if_t< Check::value > require_t
If condition is true, template is enabled.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...