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>
54 typename LLFun,
typename LLTupleArgs,
typename CovarFun,
typename CovarArgs,
56 require_t<is_all_arithmetic_scalar<CovarArgs, LLTupleArgs>>* =
nullptr>
58 CovarFun&& covariance_function,
59 CovarArgs&& covar_args,
63 [msgs, &covariance_function](
auto&&... args) {
64 return covariance_function(std::forward<
decltype(args)>(args)..., msgs);
66 std::forward<CovarArgs>(covar_args));
68 std::forward<LLFun>(ll_fun), std::forward<LLTupleArgs>(ll_args),
69 std::move(covariance), options, msgs)
75template <
bool ZeroInput =
false,
typename Output,
typename Input,
76 require_tuple_t<Output>* =
nullptr, require_tuple_t<Input>* =
nullptr,
78 std::tuple_size_v<std::decay_t<Output>> == 0>>* =
nullptr,
80 std::tuple_size_v<std::decay_t<Input>> == 0>>* =
nullptr>
91template <
bool ZeroInput =
false,
typename Output,
typename Input,
95 if constexpr (is_tuple_v<Output> && is_tuple_v<Input>) {
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 "
102 "https://github.com/stan-dev/math/issues");
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();
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();
121 }
else if constexpr (is_stan_scalar_v<output_i_t>) {
122 output_i += (0.5 * input_i.adj());
123 if constexpr (ZeroInput) {
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");
135 std::forward<Output>(output), std::forward<Input>(input));
167template <
typename LLFun,
typename LLTupleArgs,
typename CovarFun,
168 typename CovarArgs,
bool InitTheta,
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));
183 nested_rev_autodiff nested;
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);
188 [&covariance_function, &msgs](
auto&&... args) {
190 return to_var_value(covariance_function(args..., msgs));
192 return covariance_function(args..., msgs);
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);
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);
208 arena_t<Eigen::MatrixXd> LU_solve_covariance(
209 covariance.rows() * (md_est.solver_used == 3),
210 covariance.cols() * (md_est.solver_used == 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;
225 ll_fun, md_est.theta, A, options.hessian_block_size, ll_args_copy,
227 s2.deep_copy(s2_tmp);
228 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
232 * (covariance_val.diagonal() - (C.transpose() * C).diagonal())
234 ll_fun, md_est.theta, ll_args_vals, msgs))));
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;
246 options.hessian_block_size,
248 s2.deep_copy(s2_tmp);
249 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
251 }
else if (md_est.solver_used == 2) {
253 - md_est.W_r * md_est.K_root
254 * md_est.L.transpose()
255 .template triangularView<Eigen::Upper>()
257 md_est.L.template triangularView<Eigen::Lower>()
258 .solve(md_est.K_root.transpose() * md_est.W_r));
260 arena_t<Eigen::MatrixXd> C
261 = md_est.L.template triangularView<Eigen::Lower>().solve(
262 md_est.K_root.transpose());
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);
269 LU_solve_covariance = md_est.LU.solve(covariance_val);
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;
274 arena_t<Eigen::MatrixXd> A
275 = covariance_val - covariance_val * md_est.W_r * LU_solve_covariance;
277 options.hessian_block_size,
279 s2.deep_copy(s2_tmp);
280 internal::copy_compute_s2<ZeroOut>(partial_parm, ll_args_filter);
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();
292 auto covar_args_filter
296 if constexpr (ll_args_contain_var) {
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;
303 v = LU_solve_covariance * s2;
307 internal::collect_adjoints<ZeroOut>(partial_parm, ll_args_filter);
312 if constexpr (is_any_var_scalar_v<CovarArgs>) {
315 std::move(covar_args_adj));
317 if constexpr (ll_args_contain_var) {
320 std::move(partial_parm));
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.
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.
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
static void grad()
Compute the gradient for all variables starting from the end of the AD tape.
constexpr decltype(auto) apply(F &&f, Tuple &&t, PreArgs &&... pre_args)
constexpr bool is_any_var_scalar_v
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 ...