1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_ESTIMATOR_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_ESTIMATOR_HPP
16#include <unsupported/Eigen/MatrixFunctions>
67 int max_num_steps_,
bool allow_fallthrough_,
68 int max_steps_line_search_)
77template <
bool HasInitTheta>
85 hessian_block_size = hessian_block_size_;
92 Eigen::VectorXd theta_0{0};
94 template <
typename ThetaVec>
96 int hessian_block_size_,
int solver_,
97 int max_steps_line_search_,
bool allow_fallthrough_)
99 max_num_steps_, allow_fallthrough_,
100 max_steps_line_search_),
101 theta_0(
value_of(
std::forward<ThetaVec>(theta_0_))) {}
109template <
typename Options>
111 using Ops = std::decay_t<Options>;
112 if constexpr (is_tuple_v<Ops>) {
113 if constexpr (!is_eigen_v<std::tuple_element_t<0, std::decay_t<Ops>>>) {
115 sizeof(std::decay_t<Ops>*) == 0,
116 "ERROR:(laplace_marginal_lpdf) The first laplace argument is "
117 "expected to be an Eigen vector of dynamic size representing the "
120 if constexpr (!stan::is_inner_tuple_type_v<1, Ops, double>) {
122 sizeof(std::decay_t<Ops>*) == 0,
123 "ERROR:(laplace_marginal_lpdf) The second laplace argument is "
124 "expected to be a double representing the tolerance.");
126 if constexpr (!stan::is_inner_tuple_type_v<2, Ops, int>) {
128 sizeof(std::decay_t<Ops>*) == 0,
129 "ERROR:(laplace_marginal_lpdf) The third laplace argument is "
130 "expected to be an int representing the maximum number of steps for "
131 "the laplace approximation.");
133 if constexpr (!stan::is_inner_tuple_type_v<3, Ops, int>) {
135 sizeof(std::decay_t<Ops>*) == 0,
136 "ERROR:(laplace_marginal_lpdf) The fourth laplace argument is "
137 "expected to be an int representing the solver.");
139 if constexpr (!stan::is_inner_tuple_type_v<4, Ops, int>) {
141 sizeof(std::decay_t<Ops>*) == 0,
142 "ERROR:(laplace_marginal_lpdf) The fifth laplace argument is "
143 "expected to be an int representing the max steps for the laplace "
144 "approximaton's wolfe line search.");
146 constexpr bool is_fallthrough
148 5, Ops,
int> || stan::is_inner_tuple_type_v<5, Ops, bool>;
149 if constexpr (!is_fallthrough) {
151 sizeof(std::decay_t<Ops>*) == 0,
152 "ERROR:(laplace_marginal_lpdf) The sixth laplace argument is "
153 "expected to be an int representing allow fallthrough (0/1).");
157 value_of(std::get<0>(std::forward<Ops>(ops))),
160 defaults.hessian_block_size,
163 (std::get<5>(ops) > 0) ?
true :
false,
166 return std::forward<Ops>(ops);
170template <
typename ThetaVec,
typename WR,
typename L_t,
typename A_vec,
171 typename ThetaGrad,
typename LU_t,
typename KRoot>
174 double lmd{std::numeric_limits<double>::infinity()};
210 A_vec&& a_, ThetaGrad&& theta_grad_, LU_t&& LU_,
211 KRoot&& K_root_,
int solver_used_)
230template <
typename WRootMat>
232 const Eigen::SparseMatrix<double>& W,
233 const Eigen::Index block_size) {
234 int n_block = W.cols() / block_size;
235 Eigen::MatrixXd local_block(block_size, block_size);
236 Eigen::MatrixXd local_block_sqrt(block_size, block_size);
237 Eigen::MatrixXd sqrt_t_mat = Eigen::MatrixXd::Zero(block_size, block_size);
240 for (
int i = 0; i < n_block; i++) {
241 sqrt_t_mat.setZero();
243 = W.block(i * block_size, i * block_size, block_size, block_size);
244 if (!local_block.array().isFinite().any()) {
245 throw std::domain_error(
246 std::string(
"Error in block_matrix_sqrt: "
247 "NaNs detected in block diagonal starting at (")
248 + std::to_string(i) +
", " + std::to_string(i) +
")");
251 Eigen::RealSchur<Eigen::MatrixXd> schurOfA(local_block);
253 const auto& t_mat = schurOfA.matrixT();
254 const auto& u_mat = schurOfA.matrixU();
256 if ((t_mat.diagonal().array() < 0).any()) {
257 throw std::domain_error(
258 std::string(
"Error in block_matrix_sqrt: "
259 "values less than 0 detected in block diagonal's schur "
260 "decomposition starting at (")
261 + std::to_string(i) +
", " + std::to_string(i) +
")");
265 Eigen::matrix_sqrt_quasi_triangular(t_mat, sqrt_t_mat);
267 local_block_sqrt = u_mat * sqrt_t_mat * u_mat.adjoint();
268 }
catch (
const std::exception&
e) {
269 throw std::domain_error(
270 "Error in block_matrix_sqrt: "
271 "The matrix is not positive definite");
273 for (
int k = 0; k < block_size; k++) {
274 for (
int j = 0; j < block_size; j++) {
275 W_root.coeffRef(i * block_size + j, i * block_size + k)
276 = local_block_sqrt(j, k);
291template <
bool InitTheta,
typename CovarMat>
294 const CovarMat& covariance) {
295 if constexpr (InitTheta) {
297 check_finite(frame_name,
"initial guess", options.theta_0);
298 if (
unlikely(options.theta_0.size() != covariance.rows())) {
299 std::stringstream msg;
300 msg << frame_name <<
": The size of the initial theta ("
301 << options.theta_0.size()
302 <<
") vector must match the rows and columns of the covariance "
304 << covariance.rows() <<
", " << covariance.cols() <<
").";
305 throw std::domain_error(msg.str());
309 check_positive(frame_name,
"max_num_steps", options.max_num_steps);
310 check_positive(frame_name,
"hessian_block_size", options.hessian_block_size);
313 const Eigen::Index theta_size = covariance.rows();
314 if (
unlikely(theta_size % options.hessian_block_size != 0
315 || theta_size < options.hessian_block_size)) {
316 throw std::domain_error(
317 "laplace_marginal_density: Hessian block size mismatch.");
320 if (
unlikely(options.solver < 1 || options.solver > 3)) {
321 throw std::domain_error(
322 "laplace_marginal_density: solver must be 1, 2, or 3. Got: "
323 + std::to_string(options.solver));
371 template <
typename ObjFun,
typename ThetaGradFun,
typename ThetaInitializer>
372 NewtonState(
int theta_size, ObjFun&& obj_fun, ThetaGradFun&& theta_grad_f,
373 ThetaInitializer&& theta_init)
375 std::forward<ThetaInitializer>(theta_init),
376 std::forward<ThetaGradFun>(theta_grad_f)),
378 B(theta_size, theta_size),
407 template <
typename Options>
411 = std::clamp(this->
curr().alpha(), 0.0, options.line_search.max_alpha);
424template <
typename LLT,
typename B_t>
426 double max_jitter = 1
e-5) {
428 if (llt_B.info() != Eigen::Success) {
429 double jitter_try = min_jitter;
430 for (; jitter_try < max_jitter; jitter_try *= 10) {
431 B.diagonal().array() += jitter_try;
433 if (llt_B.info() == Eigen::Success) {
437 if (llt_B.info() != Eigen::Success) {
438 throw std::domain_error(
439 "laplace_marginal_density: Cholesky failed after adding jitter up to "
440 + std::to_string(jitter_try));
468 template <
typename NewtonStateT,
typename CovarMat>
490 template <
typename NewtonStateT,
typename LLFun,
typename LLTupleArgs,
493 const LLTupleArgs& ll_args,
const CovarMat& covariance,
494 int , std::ostream* msgs) {
495 const Eigen::Index theta_size = state.b.size();
500 for (Eigen::Index j = 0; j <
W_diag.size(); j++) {
501 if (
W_diag.coeff(j) < 0 || !std::isfinite(
W_diag.coeff(j))) {
502 throw std::domain_error(
503 "laplace_marginal_density: Hessian matrix is not positive "
512 = Eigen::MatrixXd::Identity(theta_size, theta_size)
518 state.b.noalias() = (
W_diag.array() * state.prev().theta().array()).matrix()
519 + state.prev().theta_grad();
520 auto L =
llt_B.matrixL();
521 auto LT =
llt_B.matrixU();
522 state.curr().a().noalias()
526 L.solve(
W_r_diag.cwiseProduct(covariance * state.b)));
534 return 2.0 *
llt_B.matrixLLT().diagonal().array().log().sum();
545 template <
typename NewtonStateT>
548 state.prev().obj() - 0.5 * log_det,
549 std::move(state).prev().
theta(),
550 Eigen::SparseMatrix<double>(
W_r_diag.asDiagonal()),
551 Eigen::MatrixXd(
llt_B.matrixL()),
552 std::move(state).prev().a(),
553 std::move(state).prev().theta_grad(),
554 Eigen::PartialPivLU<Eigen::MatrixXd>{},
555 Eigen::MatrixXd(0, 0),
575 Eigen::SparseMatrix<double>
W_r;
583 template <
typename NewtonStateT>
586 const Eigen::Index theta_size = state.b.size();
587 W_r.reserve(Eigen::VectorXi::Constant(theta_size, hessian_block_size));
588 const Eigen::Index n_block = theta_size / hessian_block_size;
589 for (Eigen::Index ii = 0; ii < n_block; ii++) {
590 for (Eigen::Index k = 0; k < hessian_block_size; k++) {
591 for (Eigen::Index j = 0; j < hessian_block_size; j++) {
592 W_r.insert(ii * hessian_block_size + j, ii * hessian_block_size + k)
597 W_r.makeCompressed();
621 template <
typename NewtonStateT,
typename LLFun,
typename LLTupleArgs,
624 const LLTupleArgs& ll_args,
const CovarMat& covariance,
625 int hessian_block_size, std::ostream* msgs) {
626 const Eigen::Index theta_size = state.b.size();
629 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
631 for (Eigen::Index j = 0; j <
W_block.rows(); j++) {
632 if (
W_block.coeff(j, j) < 0 || !std::isfinite(
W_block.coeff(j, j))) {
633 throw std::domain_error(
634 "laplace_marginal_density: Hessian matrix is not positive "
643 state.B.noalias() = Eigen::MatrixXd::Identity(theta_size, theta_size)
644 +
W_r * (covariance *
W_r);
651 =
W_block * state.prev().theta() + state.prev().theta_grad();
652 auto L =
llt_B.matrixL();
653 auto LT =
llt_B.matrixU();
654 state.curr().a().noalias()
655 = state.b -
W_r * LT.solve(L.solve(
W_r * (covariance * state.b)));
663 return 2.0 *
llt_B.matrixLLT().diagonal().array().log().sum();
674 template <
typename NewtonStateT>
677 std::move(state).prev().
theta(),
679 Eigen::MatrixXd(
llt_B.matrixL()),
680 std::move(state).prev().a(),
681 std::move(state).prev().theta_grad(),
682 Eigen::PartialPivLU<Eigen::MatrixXd>{},
683 Eigen::MatrixXd(0, 0),
710 template <
typename NewtonStateT,
typename CovarMat>
713 auto K_root_llt = covariance.template selfadjointView<Eigen::Lower>().llt();
714 if (K_root_llt.info() != Eigen::Success) {
715 throw std::domain_error(
716 "laplace_marginal_density: Cholesky of covariance failed at start");
718 K_root = std::move(K_root_llt.matrixL());
741 template <
typename NewtonStateT,
typename LLFun,
typename LLTupleArgs,
744 const LLTupleArgs& ll_args,
const CovarMat& covariance,
745 int hessian_block_size, std::ostream* msgs) {
746 const Eigen::Index theta_size = state.b.size();
750 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
753 state.B.noalias() = Eigen::MatrixXd::Identity(theta_size, theta_size)
761 =
W_full * state.prev().theta() + state.prev().theta_grad();
762 auto L =
llt_B.matrixL();
763 auto LT =
llt_B.matrixU();
764 state.curr().a().noalias()
765 =
K_root.transpose().template triangularView<Eigen::Upper>().solve(
766 LT.solve(L.solve(
K_root.transpose() * state.b)));
774 return 2.0 *
llt_B.matrixLLT().diagonal().array().log().sum();
785 template <
typename NewtonStateT>
788 std::move(state.prev().theta()),
790 Eigen::MatrixXd(
llt_B.matrixL()),
791 std::move(state.prev().a()),
792 std::move(state.prev().theta_grad()),
793 Eigen::PartialPivLU<Eigen::MatrixXd>{},
814 Eigen::PartialPivLU<Eigen::MatrixXd>
lu;
836 template <
typename NewtonStateT,
typename LLFun,
typename LLTupleArgs,
839 const LLTupleArgs& ll_args,
const CovarMat& covariance,
840 int hessian_block_size, std::ostream* msgs) {
841 const Eigen::Index theta_size = state.b.size();
845 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
848 lu.compute(Eigen::MatrixXd::Identity(theta_size, theta_size)
853 =
W_full * state.prev().theta() + state.prev().theta_grad();
854 state.curr().a().noalias()
855 = state.b -
W_full *
lu.solve(covariance * state.b);
868 return lu.matrixLU().diagonal().array().log().sum();
879 template <
typename NewtonStateT>
882 std::move(state).prev().
theta(),
884 Eigen::MatrixXd(0, 0),
885 std::move(state).prev().a(),
886 std::move(state).prev().theta_grad(),
888 Eigen::MatrixXd(0, 0),
915template <
typename SolverPolicy,
typename NewtonStateT,
typename OptionsT,
916 typename LLFunT,
typename LLTupleArgsT,
typename CovarMatT,
919 const OptionsT& options, Eigen::Index& step_iter,
920 const LLFunT& ll_fun,
const LLTupleArgsT& ll_args,
921 const CovarMatT& covariance, UpdateFun&& update_fun,
922 std::ostream* msgs) {
923 bool finish_update =
false;
924 for (; step_iter <= options.max_num_steps; step_iter++) {
925 solver.solve_step(state, ll_fun, ll_args, covariance,
926 options.hessian_block_size, msgs);
927 if (!state.final_loop) {
928 state.wolfe_info.p_ = state.curr().a() - state.prev().a();
929 state.prev_g.noalias() = -covariance * state.prev().a()
930 + covariance * state.prev().theta_grad();
931 state.wolfe_info.init_dir_ = state.prev_g.dot(state.wolfe_info.p_);
933 state.wolfe_info.flip_direction();
934 auto&& scratch = state.wolfe_info.scratch_;
935 scratch.alpha() = 1.0;
936 update_fun(scratch, state.curr(), state.prev(), scratch.eval_,
937 state.wolfe_info.p_);
938 if (scratch.alpha() <= options.line_search.min_alpha) {
939 state.wolfe_status.accept_ =
false;
940 finish_update =
true;
941 }
else if (options.line_search.max_iterations == 0) {
942 state.curr().update(scratch);
943 state.wolfe_status.accept_ =
true;
945 Eigen::VectorXd s = scratch.a() - state.prev().a();
947 = (-covariance * scratch.a() + covariance * scratch.theta_grad())
950 s, full_step_grad, state.prev_g, state.prev().alpha(),
951 state.wolfe_status.num_backtracks_, options.line_search.min_alpha,
952 options.line_search.max_alpha);
954 state.wolfe_info, update_fun, options.line_search, msgs);
960 bool objective_converged
961 = std::abs(state.curr().obj() - state.prev().obj())
963 bool search_failed = (!state.wolfe_status.accept_
964 && state.curr().obj() <= state.prev().obj());
965 finish_update = objective_converged || search_failed;
968 if (!state.final_loop && state.wolfe_status.accept_) {
970 state.final_loop =
true;
972 state.update_next_step(options);
975 return solver.build_result(state, solver.compute_log_determinant());
977 state.update_next_step(options);
983 "WARNING(laplace_marginal_density): max number of iterations: ")
984 + std::to_string(options.max_num_steps) +
" exceeded.";
986 return solver.build_result(state, solver.compute_log_determinant());
1000 std::ostream* msgs, std::string_view context,
1002 std::string_view failed_solver,
1003 std::string_view next_solver,
1004 const std::exception&
e) {
1006 std::ostringstream os;
1007 std::string msg_type = allow_fallthrough ?
"WARNING" :
"ERROR";
1008 os <<
"[" << context <<
"] " << msg_type <<
": solver fallback\n"
1009 <<
" " << std::left << std::setw(12) <<
"iteration:" << iter <<
"\n"
1010 <<
" " << std::left << std::setw(12) <<
"failed:" << failed_solver <<
"\n"
1011 <<
" " << std::left << std::setw(12) <<
"reason:" <<
e.what() <<
"\n"
1012 <<
" " << std::left << std::setw(12) <<
"action:"
1013 <<
"trying " << next_solver <<
"\n"
1014 <<
"note: this warning message will only be displayed once."
1016 if (allow_fallthrough && msgs) {
1017 (*msgs) << os.str();
1019 throw std::domain_error(std::string(
"[") + std::string(context) +
"]");
1023template <
bool InitTheta,
typename Opts>
1025 if constexpr (InitTheta) {
1027 return std::decay_t<decltype(options)>(options).theta_0;
1029 return Eigen::MatrixXd::Zero(theta_size, 1);
1051template <
typename ObjFun,
typename ThetaGradFun,
typename Covariance,
1054 Covariance&& covariance, Options&& options) {
1055 auto update_step = [&covariance, &obj_fun, &theta_grad_f](
1056 auto& proposal,
auto&& ,
auto&& prev,
1057 auto& eval_in,
auto&& p) {
1059 proposal.a() = prev.a() + eval_in.alpha() * p;
1060 proposal.theta().noalias() = covariance * proposal.a();
1061 proposal.theta_grad() = theta_grad_f(proposal.theta());
1062 eval_in.obj() = obj_fun(proposal.a(), proposal.theta());
1064 = (-covariance * proposal.a() + covariance * proposal.theta_grad())
1066 return std::isfinite(eval_in.obj()) && std::isfinite(eval_in.dir());
1067 }
catch (
const std::exception&) {
1071 auto backoff = [&options](
auto&
eval) {
1072 eval.alpha() *= options.line_search.tau;
1073 return eval.alpha() > options.line_search.min_alpha;
1076 [update_step_ = std::move(update_step), backoff_ = std::move(backoff)](
1077 auto& proposal,
auto&& curr,
auto&& prev,
auto& eval_in,
auto&& p) {
1079 eval_in, p, backoff_);
1135template <
typename LLFun,
typename LLTupleArgs,
typename CovarMat,
1139 LLFun&& ll_fun, LLTupleArgs&& ll_args, CovarMat&& covariance,
1143 const Eigen::Index theta_size = covariance.rows();
1145 auto obj_fun = [&ll_fun, &ll_args, &msgs](
const Eigen::VectorXd& a_val,
1146 auto&& theta_val) ->
double {
1147 return -0.5 * a_val.dot(theta_val)
1151 auto theta_grad_f = [&ll_fun, &ll_args, &msgs](
auto&& theta_val) {
1154 decltype(
auto) theta_init = theta_init_impl<InitTheta>(theta_size, options);
1158 std::move(obj_fun), std::move(theta_grad_f), covariance, options);
1159 Eigen::Index step_iter = 0;
1161 if (options.solver == 1) {
1162 if (options.hessian_block_size == 1) {
1165 ll_args, covariance, update_fun, msgs);
1169 ll_args, covariance, update_fun, msgs);
1172 }
catch (
const std::exception&
e) {
1173 const std::string solver_type
1174 = (options.hessian_block_size == 1) ?
"Diagonal" :
"Block";
1175 std::string failed =
"solver 1 (" + solver_type +
" Hessian-root Cholesky)";
1178 [](
auto&&... args) {
1181 options.allow_fallthrough, msgs,
"laplace_marginal_density", step_iter,
1182 std::move(failed),
"solver 2 (Covariance-root Cholesky)",
e);
1185 if (options.solver == 2 || options.allow_fallthrough) {
1187 return run_newton_loop(solver, state, options, step_iter, ll_fun, ll_args,
1188 covariance, update_fun, msgs);
1190 }
catch (
const std::exception&
e) {
1193 [](
auto&&... args) {
1196 options.allow_fallthrough, msgs,
"laplace_marginal_density", step_iter,
1197 "solver 2 (Covariance-root Cholesky)",
"solver 3 (General LU solver)",
1200 if (options.solver == 3 || options.allow_fallthrough) {
1202 return run_newton_loop(solver, state, options, step_iter, ll_fun, ll_args,
1203 covariance, update_fun, msgs);
1205 throw std::domain_error(
1206 std::string(
"You chose a solver (") + std::to_string(options.solver)
1207 +
") that is not valid. Please choose either 1, 2, or 3.");
int64_t size(const T &m)
Returns the size (number of the elements) of a matrix_cl or var_value<matrix_cl<T>>.
(Expert) Numerical traits for algorithmic differentiation variables.
WolfeStatus wolfe_line_search(Info &wolfe_info, UpdateFun &&update_fun, Options &&opt, Stream *msgs)
Strong Wolfe line search for maximization.
auto create_update_fun(ObjFun &&obj_fun, ThetaGradFun &&theta_grad_f, Covariance &&covariance, Options &&options)
Create the update function for the line search, capturing necessary references.
auto run_newton_loop(SolverPolicy &solver, NewtonStateT &state, const OptionsT &options, Eigen::Index &step_iter, const LLFunT &ll_fun, const LLTupleArgsT &ll_args, const CovarMatT &covariance, UpdateFun &&update_fun, std::ostream *msgs)
Run a Newton loop with a solver policy, updating the shared state.
double barzilai_borwein_step_size(const Eigen::VectorXd &s, const Eigen::VectorXd &g_curr, const Eigen::VectorXd &g_prev, double prev_step, int last_backtracks, double min_alpha, double max_alpha)
Curvature-aware Barzilai–Borwein (BB) step length with robust safeguards.
constexpr int laplace_default_max_num_steps
void log_solver_fallback(const bool allow_fallthrough, std::ostream *msgs, std::string_view context, Eigen::Index iter, std::string_view failed_solver, std::string_view next_solver, const std::exception &e)
Log a solver fallback event to the provided stream.
decltype(auto) theta_init_impl(Eigen::Index theta_size, Opts &&options)
constexpr double laplace_default_tolerance
constexpr int laplace_default_solver
constexpr auto tuple_to_laplace_options(Options &&ops)
constexpr int laplace_default_max_steps_line_search
constexpr int laplace_default_hessian_block_size
void validate_laplace_options(const char *frame_name, const laplace_options< InitTheta > &options, const CovarMat &covariance)
Validates the options for the Laplace approximation.
constexpr int laplace_default_allow_fallthrough
auto retry_evaluate(Update &&update, Proposal &&proposal, Curr &&curr, Prev &&prev, Eval &eval, P &&p, Backoff &&backoff)
Retry evaluation of a step until it passes a validity check.
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 llt_with_jitter(LLT &llt_B, B_t &B, double min_jitter=1e-10, double max_jitter=1e-5)
Factorize B with jittering fallback.
void block_matrix_sqrt(WRootMat &W_root, const Eigen::SparseMatrix< double > &W, const Eigen::Index block_size)
Returns the principal square root of a block diagonal matrix.
static thread_local std::once_flag fallback_warning
auto diagonal_hessian(F &&f, Theta &&theta, TupleArgs &&ll_tuple, Stream *msgs)
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.
void check_square(const char *function, const char *name, const T_y &y)
Check if the specified matrix is square.
void check_nonnegative(const char *function, const char *name, const T_y &y)
Check if y is non-negative.
static constexpr double e()
Return the base of the natural logarithm.
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 check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
void check_nonzero_size(const char *function, const char *name, const T_y &y)
Check if the specified matrix/vector is of non-zero size.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
double dot(const std::vector< double > &x, const std::vector< double > &y)
constexpr bool is_inner_tuple_type_v
Checks if the N-th element of a tuple is of the same type as CheckType.
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 ...
void solve_step(NewtonStateT &state, const LLFun &ll_fun, const LLTupleArgs &ll_args, const CovarMat &covariance, int hessian_block_size, std::ostream *msgs)
Perform one Newton step using covariance Cholesky solver.
Eigen::MatrixXd K_root
Lower Cholesky factor of covariance: Sigma = K_root * K_root^T.
Eigen::LLT< Eigen::MatrixXd > llt_B
Cholesky factorization of B = I + K_root^T * W * K_root.
Eigen::SparseMatrix< double > W_full
Full (block) Hessian matrix from likelihood.
double compute_log_determinant() const
Compute log determinant of B from Cholesky factor.
auto build_result(NewtonStateT &state, double log_det)
Build the final result structure.
CholeskyKSolver(const NewtonStateT &state, const CovarMat &covariance)
Solver Policy 2: Cholesky decomposition of K (Covariance).
Eigen::LLT< Eigen::MatrixXd > llt_B
Cholesky factorization of B = I + W_r * Sigma * W_r.
Eigen::SparseMatrix< double > W_block
Sparse block-diagonal Hessian from likelihood.
Eigen::SparseMatrix< double > W_r
Sparse square root of block Hessian.
double compute_log_determinant() const
Compute log determinant of B from Cholesky factor.
void solve_step(NewtonStateT &state, const LLFun &ll_fun, const LLTupleArgs &ll_args, const CovarMat &covariance, int hessian_block_size, std::ostream *msgs)
Perform one Newton step using block-diagonal Hessian solver.
auto build_result(NewtonStateT &state, double log_det)
Build the final result structure.
CholeskyWSolverBlock(const NewtonStateT &state, int hessian_block_size)
Solver Policy 1 (Block): Cholesky decomposition using block W.
void solve_step(NewtonStateT &state, const LLFun &ll_fun, const LLTupleArgs &ll_args, const CovarMat &covariance, int, std::ostream *msgs)
Perform one Newton step using diagonal Hessian solver.
Eigen::LLT< Eigen::MatrixXd > llt_B
Cholesky factorization of B = I + W_r * Sigma * W_r.
CholeskyWSolverDiag(const NewtonStateT &state, const CovarMat &covariance)
Eigen::VectorXd W_r_diag
Square root of diagonal Hessian: W_r[j] = sqrt(W[j])
auto build_result(NewtonStateT &state, double log_det)
Build the final result structure.
Eigen::VectorXd W_diag
Diagonal Hessian values from the likelihood.
double compute_log_determinant() const
Compute log determinant of B from Cholesky factor.
Solver Policy 1 (Diagonal): Cholesky decomposition using W.
auto build_result(NewtonStateT &state, double log_det)
Build the final result structure.
void solve_step(NewtonStateT &state, const LLFun &ll_fun, const LLTupleArgs &ll_args, const CovarMat &covariance, int hessian_block_size, std::ostream *msgs)
Perform one Newton step using LU decomposition solver.
double compute_log_determinant() const
Compute log determinant from LU factorization.
Eigen::SparseMatrix< double > W_full
Full Hessian matrix from likelihood.
Eigen::PartialPivLU< Eigen::MatrixXd > lu
LU factorization of B = I + Sigma * W.
Solver Policy 3: LU Decomposition.
auto & prev() &
Access the previous step state (mutable).
Eigen::MatrixXd B
Workspace matrix: B = I + W_r * Sigma * W_r (or similar)
WolfeStatus wolfe_status
Status of the most recent Wolfe line search.
auto & curr() &
Access the current step state (mutable).
Eigen::VectorXd b
Workspace vector: b = W * theta + grad(log_lik)
void update_next_step(const Options &options)
NewtonState(int theta_size, ObjFun &&obj_fun, ThetaGradFun &&theta_grad_f, ThetaInitializer &&theta_init)
Constructs Newton state with given dimensions and functors.
WolfeInfo wolfe_info
Wolfe line search state including current/previous steps.
const auto & curr() const &
Access the current step state (const).
Eigen::VectorXd prev_g
Previous gradient for Barzilai-Borwein step calculation.
bool final_loop
On the final loop if we found a better wolfe step, but we are going to exit, we want to make sure all...
const auto & prev() const &
Access the previous step state (const).
Holds the state for the Newton-Raphson optimization loop.
Data object used in wolfe line search.
Struct to hold the result status of the Wolfe line search.
L_t L
Solver-dependent factorization of the system matrix B.
KRoot K_root
Lower Cholesky factor of the covariance matrix.
WR W_r
Solver-dependent Hessian quantity.
A_vec a
Mode in the a parameterization, where theta = covariance * a.
ThetaGrad theta_grad
Gradient of the log-likelihood with respect to theta at the mode.
laplace_density_estimates(double lmd_, ThetaVec &&theta_, WR &&W_r_, L_t &&L_, A_vec &&a_, ThetaGrad &&theta_grad_, LU_t &&LU_, KRoot &&K_root_, int solver_used_)
Options for Wolfe line search during optimization.
laplace_options(int hessian_block_size_)
laplace_options()=default
laplace_options(ThetaVec &&theta_0_, double tolerance_, int max_num_steps_, int hessian_block_size_, int solver_, int max_steps_line_search_, bool allow_fallthrough_)
double tolerance
Iterations end when the absolute change in the optimization objective is less than this tolerance.
int solver
Which linear solver to use inside the Newton step.
laplace_options_base()=default
laplace_line_search_options line_search
laplace_options_base(int hessian_block_size_, int solver_, double tolerance_, int max_num_steps_, bool allow_fallthrough_, int max_steps_line_search_)
Options for the Laplace approximation.