1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_LIKELIHOOD_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_LIKELIHOOD_HPP
16namespace laplace_likelihood {
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);
49template <
typename F,
typename Theta,
typename Stream,
typename... Args,
51inline auto theta_grad(F&& f, Theta&& theta, Stream* msgs, Args&&... args) {
56 var f_var = f(theta_var, args..., msgs);
58 return theta_var.adj().eval();
74template <
typename F,
typename Theta,
typename Stream,
typename... Args,
76inline void ll_arg_grad(F&& f, Theta&& theta, Stream* msgs, Args&&... args) {
80 var f_var = f(theta, args..., msgs);
100template <
typename F,
typename Theta,
typename Stream,
typename... Args,
104 using Eigen::Dynamic;
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);
111 return (-hessian_v).eval();
130template <
typename F,
typename Theta,
typename Stream,
typename... Args,
133 const Eigen::Index hessian_block_size, Stream* msgs,
135 using Eigen::Dynamic;
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);
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);
148 return (-hessian_theta).eval();
151 hessian_block_size,
value_of(args)..., msgs))
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;
179 const Eigen::Index theta_size = theta.size();
180 auto theta_gradient = [&theta, &f, &msgs](
auto&&... args) {
182 Matrix<var, Dynamic, 1> theta_var = theta;
183 var f_var = f(theta_var, args..., msgs);
185 return theta_var.adj().eval();
187 if (hessian_block_size == 1) {
188 auto v = Eigen::VectorXd::Ones(theta_size);
189 Eigen::VectorXd hessian_v = Eigen::VectorXd::Zero(theta_size);
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);
197 return std::make_pair(std::move(theta_gradient), (-hessian_theta).
eval());
199 return std::make_pair(
200 std::move(theta_gradient),
220template <
typename F,
typename Theta,
typename Stream,
typename... Args,
222inline Eigen::VectorXd
third_diff(F&& f, Theta&& theta, Stream&& msgs,
225 const Eigen::Index theta_size = theta.size();
227 = std::forward<Theta>(theta);
230 for (Eigen::Index i = 0; i < theta_size; ++i) {
234 grad(ftheta_ffvar.
d_.d_.vi_);
235 return theta_var.adj().eval();
259template <
typename F,
typename Theta,
typename AMat,
typename Stream,
262 const int hessian_block_size, Stream* msgs,
264 using Eigen::Dynamic;
266 using Eigen::MatrixXd;
267 using Eigen::VectorXd;
270 const Eigen::Index theta_size = theta.size();
272 int n_blocks = theta_size / hessian_block_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) {
282 for (
int j = i; j < theta_size; j += hessian_block_size) {
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);
292 for (
int j = 0; j < theta_size; ++j) {
296 [](
auto&& f,
auto&& theta_ffvar,
auto&& msgs,
auto&&... inner_args) {
297 return f(theta_ffvar, inner_args..., msgs);
299 shallow_copy_args, f, theta_ffvar, msgs);
300 grad(target_ffvar.
d_.d_.vi_);
302 return (0.5 * theta_var.adj()).eval();
324template <
typename F,
typename V_t,
typename Theta,
typename Stream,
328 using Eigen::Dynamic;
330 using Eigen::VectorXd;
332 if constexpr (!contains_var) {
336 const Eigen::Index theta_size = theta.size();
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));
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);
348 shallow_copy_args, f, theta_fvar, msgs);
365template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
368inline auto theta_grad(F&& f, Theta&& theta, TupleArgs&& ll_tup,
369 Stream* msgs =
nullptr) {
371 [](
auto&& f,
auto&& theta,
auto&& msgs,
auto&&... args) {
373 std::forward<
decltype(theta)>(theta), msgs,
374 std::forward<
decltype(args)>(args)...);
376 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
377 std::forward<Theta>(theta), msgs);
380template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
384 Stream* msgs =
nullptr) {
386 [](
auto&& f,
auto&& theta,
auto&& msgs,
auto&&... args) {
388 std::forward<
decltype(theta)>(theta), msgs,
389 std::forward<
decltype(args)>(args)...);
391 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
392 std::forward<Theta>(theta), msgs);
395template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
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)...);
406 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
407 std::forward<Theta>(theta), msgs);
410template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
414 const Eigen::Index hessian_block_size,
415 TupleArgs&& ll_tuple, Stream* msgs) {
417 [](
auto&& f,
auto&& theta,
auto hessian_block_size,
auto* msgs,
420 std::forward<
decltype(f)>(f), std::forward<
decltype(theta)>(theta),
421 hessian_block_size, msgs, std::forward<
decltype(args)>(args)...);
423 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
424 std::forward<Theta>(theta), hessian_block_size, msgs);
438template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
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)...);
449 std::forward<TupleArgs>(ll_tup), std::forward<F>(f),
450 std::forward<Theta>(theta), msgs);
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) {
472 [](
auto&& f,
auto&& theta,
auto hessian_block_size,
auto* msgs,
475 std::forward<
decltype(f)>(f), std::forward<
decltype(theta)>(theta),
476 hessian_block_size, msgs, std::forward<
decltype(args)>(args)...);
478 std::forward<TupleArgs>(ll_tuple), std::forward<F>(f),
479 std::forward<Theta>(theta), hessian_block_size, msgs);
493template <
typename F,
typename Theta,
typename TupleArgs,
typename Stream,
496inline Eigen::VectorXd
third_diff(F&& f, Theta&& theta, TupleArgs&& ll_args,
499 [](
auto&& f,
auto&& theta,
auto&& msgs,
auto&&... args) {
501 std::forward<
decltype(theta)>(theta), msgs,
502 std::forward<
decltype(args)>(args)...);
504 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
505 std::forward<Theta>(theta), msgs);
523template <
typename F,
typename Theta,
typename AMat,
typename TupleArgs,
526inline auto compute_s2(F&& f, Theta&& theta, AMat&& A,
int hessian_block_size,
527 TupleArgs&& ll_args, Stream* msgs) {
529 [](
auto&& f,
auto&& theta,
auto&& A,
auto hessian_block_size,
auto* msgs,
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)...);
536 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
537 std::forward<Theta>(theta), std::forward<AMat>(A), hessian_block_size,
554template <
typename F,
typename V_t,
typename Theta,
typename TupleArgs,
558 TupleArgs&& ll_args, Stream* msgs) {
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)...);
566 std::forward<TupleArgs>(ll_args), std::forward<F>(f),
567 std::forward<V_t>(v), std::forward<Theta>(theta), msgs);
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.
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 ...
T value_of(const fvar< T > &v)
Return the value of the specified variable.
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.
constexpr decltype(auto) apply(F &&f, Tuple &&t, PreArgs &&... pre_args)
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.
This template class represents scalars used in forward-mode automatic differentiation,...