Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_likelihood.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_LIKELIHOOD_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_LIKELIHOOD_HPP
3
8
9namespace stan {
10namespace math {
11
16namespace laplace_likelihood {
17
18namespace internal {
29template <typename F, typename Theta, typename Stream, typename... Args,
31inline auto log_likelihood(F&& f, Theta&& theta, Stream* msgs, Args&&... args) {
32 return std::forward<F>(f)(std::forward<Theta>(theta),
33 std::forward<Args>(args)..., msgs);
34}
35
49template <typename F, typename Theta, typename Stream, typename... Args,
51inline auto theta_grad(F&& f, Theta&& theta, Stream* msgs, Args&&... args) {
52 using Eigen::Dynamic;
53 using Eigen::Matrix;
55 arena_t<Matrix<var, Dynamic, 1>> theta_var = theta;
56 var f_var = f(theta_var, args..., msgs);
57 grad(f_var.vi_);
58 return theta_var.adj().eval();
59}
60
74template <typename F, typename Theta, typename Stream, typename... Args,
76inline void ll_arg_grad(F&& f, Theta&& theta, Stream* msgs, Args&&... args) {
77 using Eigen::Dynamic;
78 using Eigen::Matrix;
80 var f_var = f(theta, args..., msgs);
81 grad(f_var.vi_);
82}
83
100template <typename F, typename Theta, typename Stream, typename... Args,
102inline auto diagonal_hessian(F&& f, Theta&& theta, Stream* msgs,
103 Args&&... args) {
104 using Eigen::Dynamic;
105 using Eigen::Matrix;
106 const Eigen::Index theta_size = theta.size();
107 auto v = Eigen::VectorXd::Ones(theta_size);
108 Eigen::VectorXd hessian_v = Eigen::VectorXd::Zero(theta_size);
109 hessian_times_vector(f, hessian_v, std::forward<Theta>(theta), std::move(v),
110 value_of(args)..., msgs);
111 return (-hessian_v).eval();
112}
113
130template <typename F, typename Theta, typename Stream, typename... Args,
132inline auto block_hessian(F&& f, Theta&& theta,
133 const Eigen::Index hessian_block_size, Stream* msgs,
134 Args&&... args) {
135 using Eigen::Dynamic;
136 using Eigen::Matrix;
137 const Eigen::Index theta_size = theta.size();
138 if (hessian_block_size == 1) {
139 auto v = Eigen::VectorXd::Ones(theta_size);
140 Eigen::VectorXd hessian_v = Eigen::VectorXd::Zero(theta_size);
141 hessian_times_vector(f, hessian_v, std::forward<Theta>(theta), std::move(v),
142 value_of(args)..., msgs);
143 Eigen::SparseMatrix<double> hessian_theta(theta_size, theta_size);
144 hessian_theta.reserve(Eigen::VectorXi::Constant(theta_size, 1));
145 for (Eigen::Index i = 0; i < theta_size; i++) {
146 hessian_theta.insert(i, i) = hessian_v(i);
147 }
148 return (-hessian_theta).eval();
149 } else {
150 return (-hessian_block_diag(f, std::forward<Theta>(theta),
151 hessian_block_size, value_of(args)..., msgs))
152 .eval();
153 }
154}
155
173template <typename F, typename Theta, typename Stream, typename... Args,
175inline auto diff(F&& f, Theta&& theta, const Eigen::Index hessian_block_size,
176 Stream* msgs, Args&&... args) {
177 using Eigen::Dynamic;
178 using Eigen::Matrix;
179 const Eigen::Index theta_size = theta.size();
180 auto theta_gradient = [&theta, &f, &msgs](auto&&... args) {
181 nested_rev_autodiff nested;
182 Matrix<var, Dynamic, 1> theta_var = theta;
183 var f_var = f(theta_var, args..., msgs);
184 grad(f_var.vi_);
185 return theta_var.adj().eval();
186 }(args...);
187 if (hessian_block_size == 1) {
188 auto v = Eigen::VectorXd::Ones(theta_size);
189 Eigen::VectorXd hessian_v = Eigen::VectorXd::Zero(theta_size);
190 hessian_times_vector(f, hessian_v, std::forward<Theta>(theta), std::move(v),
191 value_of(args)..., msgs);
192 Eigen::SparseMatrix<double> hessian_theta(theta_size, theta_size);
193 hessian_theta.reserve(Eigen::VectorXi::Constant(theta_size, 1));
194 for (Eigen::Index i = 0; i < theta_size; i++) {
195 hessian_theta.insert(i, i) = hessian_v(i);
196 }
197 return std::make_pair(std::move(theta_gradient), (-hessian_theta).eval());
198 } else {
199 return std::make_pair(
200 std::move(theta_gradient),
201 (-hessian_block_diag(f, std::forward<Theta>(theta), hessian_block_size,
202 value_of(args)..., msgs))
203 .eval());
204 }
205}
206
220template <typename F, typename Theta, typename Stream, typename... Args,
222inline Eigen::VectorXd third_diff(F&& f, Theta&& theta, Stream&& msgs,
223 Args&&... args) {
224 nested_rev_autodiff nested;
225 const Eigen::Index theta_size = theta.size();
227 = std::forward<Theta>(theta);
228 arena_t<Eigen::Matrix<fvar<fvar<var>>, Eigen::Dynamic, 1>> theta_ffvar(
229 theta_size);
230 for (Eigen::Index i = 0; i < theta_size; ++i) {
231 theta_ffvar(i) = fvar<fvar<var>>(fvar<var>(theta_var(i), 1.0), 1.0);
232 }
233 fvar<fvar<var>> ftheta_ffvar = f(theta_ffvar, args..., msgs);
234 grad(ftheta_ffvar.d_.d_.vi_);
235 return theta_var.adj().eval();
236}
237
259template <typename F, typename Theta, typename AMat, typename Stream,
260 typename... Args, require_eigen_vector_t<Theta>* = nullptr>
261inline auto compute_s2(F&& f, Theta&& theta, AMat&& A,
262 const int hessian_block_size, Stream* msgs,
263 Args&&... args) {
264 using Eigen::Dynamic;
265 using Eigen::Matrix;
266 using Eigen::MatrixXd;
267 using Eigen::VectorXd;
268
269 nested_rev_autodiff nested;
270 const Eigen::Index theta_size = theta.size();
271 arena_t<Matrix<var, Dynamic, 1>> theta_var = std::forward<Theta>(theta);
272 int n_blocks = theta_size / hessian_block_size;
273 arena_t<VectorXd> v(theta_size);
274 arena_t<VectorXd> w(theta_size);
275 Matrix<fvar<fvar<var>>, Dynamic, 1> theta_ffvar(theta_size);
276 auto shallow_copy_args
277 = stan::math::internal::shallow_copy_vargs<fvar<fvar<var>>>(
278 std::forward_as_tuple(args...));
279 for (Eigen::Index i = 0; i < hessian_block_size; ++i) {
280 nested_rev_autodiff nested;
281 v.setZero();
282 for (int j = i; j < theta_size; j += hessian_block_size) {
283 v(j) = 1;
284 }
285 w.setZero();
286 for (int j = 0; j < n_blocks; ++j) {
287 for (int k = 0; k < hessian_block_size; ++k) {
288 w(k + j * hessian_block_size)
289 = A(k + j * hessian_block_size, i + j * hessian_block_size);
290 }
291 }
292 for (int j = 0; j < theta_size; ++j) {
293 theta_ffvar(j) = fvar<fvar<var>>(fvar<var>(theta_var(j), v(j)), w(j));
294 }
295 fvar<fvar<var>> target_ffvar = stan::math::apply(
296 [](auto&& f, auto&& theta_ffvar, auto&& msgs, auto&&... inner_args) {
297 return f(theta_ffvar, inner_args..., msgs);
298 },
299 shallow_copy_args, f, theta_ffvar, msgs);
300 grad(target_ffvar.d_.d_.vi_);
301 }
302 return (0.5 * theta_var.adj()).eval();
303}
304
324template <typename F, typename V_t, typename Theta, typename Stream,
325 typename... Args, require_eigen_vector_t<Theta>* = nullptr>
326inline auto diff_eta_implicit(F&& f, V_t&& v, Theta&& theta, Stream* msgs,
327 Args&&... args) {
328 using Eigen::Dynamic;
329 using Eigen::Matrix;
330 using Eigen::VectorXd;
331 constexpr bool contains_var = is_any_var_scalar<Args...>::value;
332 if constexpr (!contains_var) {
333 return;
334 }
335 nested_rev_autodiff nested;
336 const Eigen::Index theta_size = theta.size();
337 arena_t<Matrix<var, Dynamic, 1>> theta_var = std::forward<Theta>(theta);
338 Matrix<fvar<var>, Dynamic, 1> theta_fvar(theta_size);
339 for (Eigen::Index i = 0; i < theta_size; i++) {
340 theta_fvar(i) = fvar<var>(theta_var(i), v(i));
341 }
342 auto shallow_copy_args = stan::math::internal::shallow_copy_vargs<fvar<var>>(
343 std::forward_as_tuple(args...));
345 [](auto&& f, auto&& theta_fvar, auto&& msgs, auto&&... inner_args) {
346 return f(theta_fvar, inner_args..., msgs);
347 },
348 shallow_copy_args, f, theta_fvar, msgs);
349 grad(f_sum.d_.vi_);
350}
351
352} // namespace internal
353
365template <typename F, typename Theta, typename TupleArgs, typename Stream,
368inline auto theta_grad(F&& f, Theta&& theta, TupleArgs&& ll_tup,
369 Stream* msgs = nullptr) {
370 return apply(
371 [](auto&& f, auto&& theta, auto&& msgs, auto&&... args) {
372 return internal::theta_grad(std::forward<decltype(f)>(f),
373 std::forward<decltype(theta)>(theta), msgs,
374 std::forward<decltype(args)>(args)...);
375 },
376 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
377 std::forward<Theta>(theta), msgs);
378}
379
380template <typename F, typename Theta, typename TupleArgs, typename Stream,
383inline auto ll_arg_grad(F&& f, Theta&& theta, TupleArgs&& ll_tup,
384 Stream* msgs = nullptr) {
385 return apply(
386 [](auto&& f, auto&& theta, auto&& msgs, auto&&... args) {
387 return internal::ll_arg_grad(std::forward<decltype(f)>(f),
388 std::forward<decltype(theta)>(theta), msgs,
389 std::forward<decltype(args)>(args)...);
390 },
391 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
392 std::forward<Theta>(theta), msgs);
393}
394
395template <typename F, typename Theta, typename TupleArgs, typename Stream,
398inline auto diagonal_hessian(F&& f, Theta&& theta, TupleArgs&& ll_tuple,
399 Stream* msgs) {
400 return apply(
401 [](auto&& f, auto&& theta, auto* msgs, auto&&... args) {
403 std::forward<decltype(f)>(f), std::forward<decltype(theta)>(theta),
404 msgs, std::forward<decltype(args)>(args)...);
405 },
406 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
407 std::forward<Theta>(theta), msgs);
408}
409
410template <typename F, typename Theta, typename TupleArgs, typename Stream,
413inline auto block_hessian(F&& f, Theta&& theta,
414 const Eigen::Index hessian_block_size,
415 TupleArgs&& ll_tuple, Stream* msgs) {
416 return apply(
417 [](auto&& f, auto&& theta, auto hessian_block_size, auto* msgs,
418 auto&&... args) {
420 std::forward<decltype(f)>(f), std::forward<decltype(theta)>(theta),
421 hessian_block_size, msgs, std::forward<decltype(args)>(args)...);
422 },
423 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
424 std::forward<Theta>(theta), hessian_block_size, msgs);
425}
426
438template <typename F, typename Theta, typename TupleArgs, typename Stream,
441inline auto log_likelihood(F&& f, Theta&& theta, TupleArgs&& ll_tup,
442 Stream* msgs) {
443 return apply(
444 [](auto&& f, auto&& theta, auto&& msgs, auto&&... args) {
446 std::forward<decltype(f)>(f), std::forward<decltype(theta)>(theta),
447 msgs, std::forward<decltype(args)>(args)...);
448 },
449 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
450 std::forward<Theta>(theta), msgs);
451}
452
466template <typename F, typename Theta, typename TupleArgs, typename Stream,
469inline auto diff(F&& f, Theta&& theta, const Eigen::Index hessian_block_size,
470 TupleArgs&& ll_tuple, Stream* msgs) {
471 return apply(
472 [](auto&& f, auto&& theta, auto hessian_block_size, auto* msgs,
473 auto&&... args) {
474 return internal::diff(
475 std::forward<decltype(f)>(f), std::forward<decltype(theta)>(theta),
476 hessian_block_size, msgs, std::forward<decltype(args)>(args)...);
477 },
478 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
479 std::forward<Theta>(theta), hessian_block_size, msgs);
480}
481
493template <typename F, typename Theta, typename TupleArgs, typename Stream,
496inline Eigen::VectorXd third_diff(F&& f, Theta&& theta, TupleArgs&& ll_args,
497 Stream* msgs) {
498 return apply(
499 [](auto&& f, auto&& theta, auto&& msgs, auto&&... args) {
500 return internal::third_diff(std::forward<decltype(f)>(f),
501 std::forward<decltype(theta)>(theta), msgs,
502 std::forward<decltype(args)>(args)...);
503 },
504 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
505 std::forward<Theta>(theta), msgs);
506}
507
523template <typename F, typename Theta, typename AMat, typename TupleArgs,
524 typename Stream, require_eigen_vector_t<Theta>* = nullptr,
526inline auto compute_s2(F&& f, Theta&& theta, AMat&& A, int hessian_block_size,
527 TupleArgs&& ll_args, Stream* msgs) {
528 return apply(
529 [](auto&& f, auto&& theta, auto&& A, auto hessian_block_size, auto* msgs,
530 auto&&... args) {
532 std::forward<decltype(f)>(f), std::forward<decltype(theta)>(theta),
533 std::forward<decltype(A)>(A), hessian_block_size, msgs,
534 std::forward<decltype(args)>(args)...);
535 },
536 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
537 std::forward<Theta>(theta), std::forward<AMat>(A), hessian_block_size,
538 msgs);
539}
540
554template <typename F, typename V_t, typename Theta, typename TupleArgs,
555 typename Stream, require_tuple_t<TupleArgs>* = nullptr,
557inline auto diff_eta_implicit(F&& f, V_t&& v, Theta&& theta,
558 TupleArgs&& ll_args, Stream* msgs) {
559 return apply(
560 [](auto&& f, auto&& v, auto&& theta, auto&& msgs, auto&&... args) {
562 std::forward<decltype(f)>(f), std::forward<decltype(v)>(v),
563 std::forward<decltype(theta)>(theta), msgs,
564 std::forward<decltype(args)>(args)...);
565 },
566 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
567 std::forward<V_t>(v), std::forward<Theta>(theta), msgs);
568}
569
570} // namespace laplace_likelihood
571
572} // namespace math
573} // namespace stan
574
575#endif
A class following the RAII idiom to start and recover nested autodiff scopes.
require_t< is_eigen_vector< std::decay_t< T > > > require_eigen_vector_t
Require type satisfies is_eigen_vector.
require_t< container_type_check_base< is_eigen_vector, value_type_t, TypeCheck, Check... > > require_eigen_vector_vt
Require type satisfies is_eigen_vector.
require_t< is_tuple< std::decay_t< T > > > require_tuple_t
Require type satisfies is_tuple.
Definition is_tuple.hpp:33
auto block_hessian(F &&f, Theta &&theta, const Eigen::Index hessian_block_size, Stream *msgs, Args &&... args)
Computes negative block diagonal Hessian of f wrttheta and args...
Eigen::VectorXd third_diff(F &&f, Theta &&theta, Stream &&msgs, Args &&... args)
Compute third order derivative of f wrt theta and args...
auto diff_eta_implicit(F &&f, V_t &&v, Theta &&theta, Stream *msgs, Args &&... args)
Compute second order gradient of f wrt theta and args...
auto log_likelihood(F &&f, Theta &&theta, Stream *msgs, Args &&... args)
auto diff(F &&f, Theta &&theta, const Eigen::Index hessian_block_size, Stream *msgs, Args &&... args)
Computes theta gradient and negative block diagonal Hessian of f wrt theta and args....
auto diagonal_hessian(F &&f, Theta &&theta, Stream *msgs, Args &&... args)
Computes negative diagonal Hessian of f wrttheta and args...
void ll_arg_grad(F &&f, Theta &&theta, Stream *msgs, Args &&... args)
Computes likelihood argument gradient of f
auto compute_s2(F &&f, Theta &&theta, AMat &&A, const int hessian_block_size, Stream *msgs, Args &&... args)
The derivative of the log likelihood wrt theta evaluated at the mode.
auto theta_grad(F &&f, Theta &&theta, Stream *msgs, Args &&... args)
Computes theta gradient f wrt theta and args...
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(F &&f, Theta &&theta, const Eigen::Index hessian_block_size, TupleArgs &&ll_tuple, Stream *msgs)
A wrapper that accepts a tuple as arguments.
auto diagonal_hessian(F &&f, Theta &&theta, TupleArgs &&ll_tuple, Stream *msgs)
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 log_likelihood(F &&f, Theta &&theta, TupleArgs &&ll_tup, Stream *msgs)
A wrapper that accepts a tuple as arguments.
auto block_hessian(F &&f, Theta &&theta, const Eigen::Index hessian_block_size, TupleArgs &&ll_tuple, Stream *msgs)
auto theta_grad(F &&f, Theta &&theta, TupleArgs &&ll_tup, Stream *msgs=nullptr)
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.
Eigen::SparseMatrix< double > hessian_block_diag(F &&f, const Eigen::VectorXd &x, const Eigen::Index hessian_block_size, Args &&... args)
Returns a block diagonal Hessian by computing the relevant directional derivatives and storing them i...
T eval(T &&arg)
Inputs which have a plain_type equal to the own time are forwarded unmodified (for Eigen expressions ...
Definition eval.hpp:20
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition value_of.hpp:18
void hessian_times_vector(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &Hv)
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
typename internal::arena_type_impl< std::decay_t< T > >::type arena_t
Determines a type that can be used in place of T that does any dynamic allocations on the AD stack.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
Scalar d_
The tangent (derivative) of this variable.
Definition fvar.hpp:61
This template class represents scalars used in forward-mode automatic differentiation,...
Definition fvar.hpp:40