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
148template <typename Output, typename Input>
149inline void reverse_pass_collect_adjoints(var ret, Output&& output,
150 Input&& input) {
151 if constexpr (is_tuple_v<Output>) {
153 [ret](auto&& inner_arg, auto&& inner_input) mutable {
155 ret, std::forward<decltype(inner_arg)>(inner_arg),
156 std::forward<decltype(inner_input)>(inner_input));
157 },
158 std::forward<Output>(output), std::forward<Input>(input));
159 } else if constexpr (is_std_vector_containing_tuple_v<Output>) {
160 for (std::size_t i = 0; i < output.size(); ++i) {
161 reverse_pass_collect_adjoints(ret, output[i], input[i]);
162 }
163 } else {
165 [vi = ret.vi_, arg_arena = to_arena(std::forward<Output>(output)),
166 input_arena = to_arena(std::forward<Input>(input))]() mutable {
167 collect_adjoints(arg_arena, vi, input_arena);
168 });
169 }
170}
171
172} // namespace internal
173
201template <typename LLFun, typename LLTupleArgs, typename CovarFun,
202 typename CovarArgs, bool InitTheta,
204inline auto laplace_marginal_density(LLFun&& ll_fun, LLTupleArgs&& ll_args,
205 CovarFun&& covariance_function,
206 CovarArgs&& covar_args,
207 const laplace_options<InitTheta>& options,
208 std::ostream* msgs) {
209 auto covar_args_refs = to_ref(std::forward<CovarArgs>(covar_args));
210 auto ll_args_refs = to_ref(std::forward<LLTupleArgs>(ll_args));
211 // Solver 1, 2, 3
212 constexpr bool ll_args_contain_var = is_any_var_scalar<LLTupleArgs>::value;
213 auto partial_parm = internal::make_zeroed_arena(ll_args_refs);
214 auto covar_args_adj = internal::make_zeroed_arena(covar_args_refs);
215 double lmd = 0.0;
216 {
217 nested_rev_autodiff nested;
218 // Make one hard copy here
219 auto ll_args_copy = internal::deep_copy_vargs<var>(ll_args_refs);
220 auto covar_args_copy = internal::deep_copy_vargs<var>(covar_args_refs);
221 auto covariance = stan::math::apply(
222 [&covariance_function, &msgs](auto&&... args) {
223 if constexpr (is_any_var_scalar_v<decltype(args)...>) {
224 return to_var_value(covariance_function(args..., msgs));
225 } else {
226 return covariance_function(args..., msgs);
227 }
228 },
229 covar_args_copy);
230 decltype(auto) covariance_val = value_of(covariance);
231 decltype(auto) ll_args_vals = value_of(ll_args_copy);
233 ll_fun, ll_args_vals, covariance_val, options, msgs);
234 auto ll_args_filter = internal::filter_var_scalar_types(ll_args_copy);
235 // tuple of references to var types
236 // Solver 1, 2
237 const bool solver_1_or_2
238 = md_est.solver_used == 1 || md_est.solver_used == 2;
239 arena_t<Eigen::MatrixXd> R(md_est.theta.size() * solver_1_or_2,
240 md_est.theta.size() * solver_1_or_2);
241 // Solver 3
242 arena_t<Eigen::MatrixXd> LU_solve_covariance(
243 covariance.rows() * (md_est.solver_used == 3),
244 covariance.cols() * (md_est.solver_used == 3));
245 // Solver 1, 2, 3
246 arena_t<Eigen::VectorXd> s2(md_est.theta.size());
248 if (md_est.solver_used == 1) {
249 if (options.hessian_block_size == 1) {
250 arena_t<Eigen::MatrixXd> tmp = md_est.W_r.toDense();
251 md_est.L.template triangularView<Eigen::Lower>().solveInPlace(tmp);
252 R.noalias() = tmp.transpose() * tmp;
253 arena_t<Eigen::MatrixXd> C
254 = md_est.L.template triangularView<Eigen::Lower>().solve(
255 md_est.W_r * covariance_val);
256 if constexpr (ll_args_contain_var) {
257 arena_t<Eigen::MatrixXd> A = covariance_val - C.transpose() * C;
258 auto s2_tmp = laplace_likelihood::compute_s2(
259 ll_fun, md_est.theta, A, options.hessian_block_size, ll_args_copy,
260 msgs);
261 s2.deep_copy(s2_tmp);
262 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
263 } else {
264 s2.deep_copy(
265 (0.5
266 * (covariance_val.diagonal() - (C.transpose() * C).diagonal())
267 .cwiseProduct(laplace_likelihood::third_diff(
268 ll_fun, md_est.theta, ll_args_vals, msgs))));
269 }
270
271 } else {
272 arena_t<Eigen::MatrixXd> tmp = md_est.W_r.toDense();
273 md_est.L.template triangularView<Eigen::Lower>().solveInPlace(tmp);
274 R.noalias() = tmp.transpose() * tmp;
275 arena_t<Eigen::MatrixXd> C
276 = md_est.L.template triangularView<Eigen::Lower>().solve(
277 md_est.W_r * covariance_val);
278 arena_t<Eigen::MatrixXd> A = covariance_val - C.transpose() * C;
279 auto s2_tmp = laplace_likelihood::compute_s2(ll_fun, md_est.theta, A,
280 options.hessian_block_size,
281 ll_args_copy, msgs);
282 s2.deep_copy(s2_tmp);
283 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
284 }
285 } else if (md_est.solver_used == 2) {
286 R = md_est.W_r
287 - md_est.W_r * md_est.K_root
288 * md_est.L.transpose()
289 .template triangularView<Eigen::Upper>()
290 .solve(
291 md_est.L.template triangularView<Eigen::Lower>()
292 .solve(md_est.K_root.transpose() * md_est.W_r));
293
294 arena_t<Eigen::MatrixXd> C
295 = md_est.L.template triangularView<Eigen::Lower>().solve(
296 md_est.K_root.transpose());
297 auto s2_tmp = laplace_likelihood::compute_s2(
298 ll_fun, md_est.theta, (C.transpose() * C).eval(),
299 options.hessian_block_size, ll_args_copy, msgs);
300 s2.deep_copy(s2_tmp);
301 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
302 } else { // options.solver with LU decomposition
303 LU_solve_covariance = md_est.LU.solve(covariance_val);
304 auto I_minus_BinvKW
305 = Eigen::MatrixXd::Identity(md_est.W_r.rows(), md_est.W_r.cols())
306 - LU_solve_covariance * md_est.W_r;
307 R = md_est.W_r * I_minus_BinvKW; // == W - W B^{-1} K W
308 arena_t<Eigen::MatrixXd> A
309 = covariance_val - covariance_val * md_est.W_r * LU_solve_covariance;
310 auto s2_tmp = laplace_likelihood::compute_s2(ll_fun, md_est.theta, A,
311 options.hessian_block_size,
312 ll_args_copy, msgs);
313 s2.deep_copy(s2_tmp);
314 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
315 }
316 if constexpr (is_any_var_scalar_v<scalar_type_t<CovarArgs>>) {
317 arena_t<Eigen::MatrixXd> K_adj_arena
318 = 0.5 * md_est.a * md_est.a.transpose() - 0.5 * R
319 + s2 * md_est.theta_grad.transpose()
320 - (R * (covariance.val() * s2)) * md_est.theta_grad.transpose();
322 0.0, [covariance, K_adj_arena](auto&& vi) mutable {
323 covariance.adj().array() += vi.adj() * K_adj_arena.array();
324 });
325 grad(Z.vi_);
326 auto covar_args_filter
327 = internal::filter_var_scalar_types(covar_args_copy);
328 internal::collect_adjoints(covar_args_adj, covar_args_filter);
329 }
330 if constexpr (ll_args_contain_var) {
331 laplace_likelihood::ll_arg_grad(ll_fun, md_est.theta, ll_args_copy, msgs);
332 internal::collect_adjoints<ZeroOut>(partial_parm, ll_args_filter);
333 arena_t<Eigen::VectorXd> v;
334 if (md_est.solver_used == 1 || md_est.solver_used == 2) {
335 v = covariance_val * s2 - covariance_val * R * covariance_val * s2;
336 } else {
337 v = LU_solve_covariance * s2;
338 }
339 laplace_likelihood::diff_eta_implicit(ll_fun, v, md_est.theta,
340 ll_args_copy, msgs);
341 internal::collect_adjoints<ZeroOut>(partial_parm, ll_args_filter);
342 }
343 lmd = md_est.lmd;
344 }
345 var ret(lmd);
346 if constexpr (is_any_var_scalar_v<CovarArgs>) {
347 auto covar_args_filter = internal::filter_var_scalar_types(covar_args_refs);
348 internal::reverse_pass_collect_adjoints(ret, covar_args_filter,
349 std::move(covar_args_adj));
350 }
351 if constexpr (ll_args_contain_var) {
352 auto ll_args_filter = internal::filter_var_scalar_types(ll_args_refs);
353 internal::reverse_pass_collect_adjoints(ret, ll_args_filter,
354 std::move(partial_parm));
355 }
356 return ret;
357}
358
359} // namespace math
360} // namespace stan
361
362#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...
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
constexpr void for_each(F &&f, const std::tuple<> &)
Apply a function to each element of a tuple.
Definition for_each.hpp:80
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,...
arena_t< T > to_arena(const T &a)
Converts given argument into a type that either has any dynamic allocation on AD stack or schedules i...
Definition to_arena.hpp:25
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 ...