Automatic Differentiation
 
Loading...
Searching...
No Matches
laplace_marginal_density_estimator.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_ESTIMATOR_HPP
2#define STAN_MATH_MIX_FUNCTOR_LAPLACE_MARGINAL_DENSITY_ESTIMATOR_HPP
14#include <unsupported/Eigen/MatrixFunctions>
15#include <cmath>
16#include <mutex>
17
25namespace stan {
26namespace math {
27
32 /* Size of the blocks in block diagonal hessian*/
33 int hessian_block_size{1}; // 0
49 int solver{1}; // 1
57 double tolerance{1.49012e-08}; // 2
58 /* Maximum number of steps*/
59 int max_num_steps{500}; // 3
60 int allow_fallthrough{true}; // 4
63 laplace_options_base(int hessian_block_size_, int solver_, double tolerance_,
64 int max_num_steps_, bool allow_fallthrough_,
65 int max_steps_line_search_)
66 : hessian_block_size(hessian_block_size_),
67 solver(solver_),
68 tolerance(tolerance_),
69 max_num_steps(max_num_steps_),
70 allow_fallthrough(allow_fallthrough_),
71 line_search(max_steps_line_search_) {}
72};
73
74template <bool HasInitTheta>
76
77template <>
78struct laplace_options<false> : public laplace_options_base {};
79
80template <>
82 /* Value for user supplied initial theta */
83 Eigen::VectorXd theta_0{0}; // 6
84
85 template <typename ThetaVec>
86 laplace_options(ThetaVec&& theta_0_, double tolerance_, int max_num_steps_,
87 int hessian_block_size_, int solver_,
88 int max_steps_line_search_, bool allow_fallthrough_)
89 : laplace_options_base(hessian_block_size_, solver_, tolerance_,
90 max_num_steps_, allow_fallthrough_,
91 max_steps_line_search_),
92 theta_0(std::forward<ThetaVec>(theta_0_)) {}
93};
94
97
103inline auto generate_laplace_options(int theta_0_size) {
104 auto ops = laplace_options_default{};
105 return std::make_tuple(
106 Eigen::VectorXd::Zero(theta_0_size).eval(), // 0 -> 6
107 ops.tolerance, ops.max_num_steps, ops.hessian_block_size, ops.solver,
108 ops.line_search.max_iterations, static_cast<int>(ops.allow_fallthrough));
109}
110
111namespace internal {
112
113template <typename Options>
114inline constexpr auto tuple_to_laplace_options(Options&& ops) {
115 using Ops = std::decay_t<Options>;
116 if constexpr (is_tuple_v<Ops>) {
117 if constexpr (!is_eigen_v<std::tuple_element_t<0, std::decay_t<Ops>>>) {
118 static_assert(
119 sizeof(std::decay_t<Ops>*) == 0,
120 "ERROR:(laplace_marginal_lpdf) The first laplace argument is "
121 "expected to be an Eigen vector of dynamic size representing the "
122 "initial theta_0.");
123 }
124 if constexpr (!stan::is_inner_tuple_type_v<1, Ops, double>) {
125 static_assert(
126 sizeof(std::decay_t<Ops>*) == 0,
127 "ERROR:(laplace_marginal_lpdf) The second laplace argument is "
128 "expected to be a double representing the tolerance.");
129 }
130 if constexpr (!stan::is_inner_tuple_type_v<2, Ops, int>) {
131 static_assert(
132 sizeof(std::decay_t<Ops>*) == 0,
133 "ERROR:(laplace_marginal_lpdf) The third laplace argument is "
134 "expected to be an int representing the maximum number of steps for "
135 "the laplace approximation.");
136 }
137 if constexpr (!stan::is_inner_tuple_type_v<3, Ops, int>) {
138 static_assert(
139 sizeof(std::decay_t<Ops>*) == 0,
140 "ERROR:(laplace_marginal_lpdf) The fifth laplace argument is "
141 "expected to be an int representing the hessian block size.");
142 }
143 if constexpr (!stan::is_inner_tuple_type_v<4, Ops, int>) {
144 static_assert(
145 sizeof(std::decay_t<Ops>*) == 0,
146 "ERROR:(laplace_marginal_lpdf) The fourth laplace argument is "
147 "expected to be an int representing the solver.");
148 }
149 if constexpr (!stan::is_inner_tuple_type_v<5, Ops, int>) {
150 static_assert(
151 sizeof(std::decay_t<Ops>*) == 0,
152 "ERROR:(laplace_marginal_lpdf) The sixth laplace argument is "
153 "expected to be an int representing the max steps for the laplace "
154 "approximaton's wolfe line search.");
155 }
156 constexpr bool is_fallthrough
158 6, Ops, int> || stan::is_inner_tuple_type_v<6, Ops, bool>;
159 if constexpr (!is_fallthrough) {
160 static_assert(
161 sizeof(std::decay_t<Ops>*) == 0,
162 "ERROR:(laplace_marginal_lpdf) The seventh laplace argument is "
163 "expected to be an int representing allow fallthrough (0/1).");
164 }
166 value_of(std::get<0>(std::forward<Ops>(ops))),
167 std::get<1>(ops),
168 std::get<2>(ops),
169 std::get<3>(ops),
170 std::get<4>(ops),
171 std::get<5>(ops),
172 (std::get<6>(ops) > 0) ? true : false,
173 };
174 } else {
175 return std::forward<Ops>(ops);
176 }
177}
178
179template <typename ThetaVec, typename WR, typename L_t, typename A_vec,
180 typename ThetaGrad, typename LU_t, typename KRoot>
182 /* log marginal density */
183 double lmd{std::numeric_limits<double>::infinity()};
184 /* ThetaVec at the mode */
185 ThetaVec theta;
193 WR W_r;
200 L_t L;
205 A_vec a;
207 ThetaGrad theta_grad;
208 /* LU matrix from solver 3 */
209 LU_t LU;
216 KRoot K_root;
218 laplace_density_estimates(double lmd_, ThetaVec&& theta_, WR&& W_r_, L_t&& L_,
219 A_vec&& a_, ThetaGrad&& theta_grad_, LU_t&& LU_,
220 KRoot&& K_root_, int solver_used_)
221 : lmd(lmd_),
222 theta(std::move(theta_)),
223 W_r(std::move(W_r_)),
224 L(std::move(L_)),
225 a(std::move(a_)),
226 theta_grad(std::move(theta_grad_)),
227 LU(std::move(LU_)),
228 K_root(std::move(K_root_)),
229 solver_used(solver_used_) {}
230};
231
239template <typename WRootMat>
240inline void block_matrix_sqrt(WRootMat& W_root,
241 const Eigen::SparseMatrix<double>& W,
242 const Eigen::Index block_size) {
243 int n_block = W.cols() / block_size;
244 Eigen::MatrixXd local_block(block_size, block_size);
245 Eigen::MatrixXd local_block_sqrt(block_size, block_size);
246 Eigen::MatrixXd sqrt_t_mat = Eigen::MatrixXd::Zero(block_size, block_size);
247 // No block operation available for sparse matrices, so we have to loop
248 // See https://eigen.tuxfamily.org/dox/group__TutorialSparse.html#title7
249 for (int i = 0; i < n_block; i++) {
250 sqrt_t_mat.setZero();
251 local_block
252 = W.block(i * block_size, i * block_size, block_size, block_size);
253 if (!local_block.array().isFinite().any()) {
254 throw std::domain_error(
255 std::string("Error in block_matrix_sqrt: "
256 "NaNs detected in block diagonal starting at (")
257 + std::to_string(i) + ", " + std::to_string(i) + ")");
258 }
259 // Issue here, sqrt is done over T of the complex schur
260 Eigen::RealSchur<Eigen::MatrixXd> schurOfA(local_block);
261 // Compute Schur decomposition of arg
262 const auto& t_mat = schurOfA.matrixT();
263 const auto& u_mat = schurOfA.matrixU();
264 // Check if diagonal of schur is not positive
265 if ((t_mat.diagonal().array() < 0).any()) {
266 throw std::domain_error(
267 std::string("Error in block_matrix_sqrt: "
268 "values less than 0 detected in block diagonal's schur "
269 "decomposition starting at (")
270 + std::to_string(i) + ", " + std::to_string(i) + ")");
271 }
272 try {
273 // Compute square root of T
274 Eigen::matrix_sqrt_quasi_triangular(t_mat, sqrt_t_mat);
275 // Compute square root of arg
276 local_block_sqrt = u_mat * sqrt_t_mat * u_mat.adjoint();
277 } catch (const std::exception& e) {
278 throw std::domain_error(
279 "Error in block_matrix_sqrt: "
280 "The matrix is not positive definite");
281 }
282 for (int k = 0; k < block_size; k++) {
283 for (int j = 0; j < block_size; j++) {
284 W_root.coeffRef(i * block_size + j, i * block_size + k)
285 = local_block_sqrt(j, k);
286 }
287 }
288 }
289}
290
300template <bool InitTheta, typename CovarMat>
301inline void validate_laplace_options(const char* frame_name,
302 const laplace_options<InitTheta>& options,
303 const CovarMat& covariance) {
304 if constexpr (InitTheta) {
305 check_nonzero_size(frame_name, "initial guess", options.theta_0);
306 check_finite(frame_name, "initial guess", options.theta_0);
307 if (unlikely(options.theta_0.size() != covariance.rows())) {
308 std::stringstream msg;
309 msg << frame_name << ": The size of the initial theta ("
310 << options.theta_0.size()
311 << ") vector must match the rows and columns of the covariance "
312 "matrix ("
313 << covariance.rows() << ", " << covariance.cols() << ").";
314 throw std::domain_error(msg.str());
315 }
316 }
317 check_nonnegative(frame_name, "tolerance", options.tolerance);
318 check_positive(frame_name, "max_num_steps", options.max_num_steps);
319 check_positive(frame_name, "hessian_block_size", options.hessian_block_size);
320 check_square(frame_name, "covariance", covariance);
321
322 const Eigen::Index theta_size = covariance.rows();
323 if (unlikely(theta_size % options.hessian_block_size != 0
324 || theta_size < options.hessian_block_size)) {
325 throw std::domain_error(
326 "laplace_marginal_density: Hessian block size mismatch.");
327 }
328
329 if (unlikely(options.solver < 1 || options.solver > 3)) {
330 throw std::domain_error(
331 "laplace_marginal_density: solver must be 1, 2, or 3. Got: "
332 + std::to_string(options.solver));
333 }
334}
335
350
353
355 Eigen::VectorXd b;
356
358 Eigen::MatrixXd B;
359
361 Eigen::VectorXd prev_g;
369 bool final_loop = false;
370
380 template <typename ObjFun, typename ThetaGradFun, typename ThetaInitializer>
381 NewtonState(int theta_size, ObjFun&& obj_fun, ThetaGradFun&& theta_grad_f,
382 ThetaInitializer&& theta_init)
383 : wolfe_info(std::forward<ObjFun>(obj_fun), theta_size,
384 std::forward<ThetaInitializer>(theta_init),
385 std::forward<ThetaGradFun>(theta_grad_f)),
386 b(theta_size),
387 B(theta_size, theta_size),
388 prev_g(theta_size) {
389 wolfe_status.num_backtracks_ = -1; // Safe initial value for BB step
390 }
391
396 auto& curr() & { return wolfe_info.curr_; }
397
402 const auto& curr() const& { return wolfe_info.curr_; }
403 auto&& curr() && { return std::move(wolfe_info).curr(); }
408 auto& prev() & { return wolfe_info.prev_; }
409
414 const auto& prev() const& { return wolfe_info.prev_; }
415 auto&& prev() && { return std::move(wolfe_info).prev(); }
416 template <typename Options>
417 inline void update_next_step(const Options& options) {
418 this->prev().update(this->curr());
419 this->curr().alpha()
420 = std::clamp(this->curr().alpha(), 0.0, options.line_search.max_alpha);
421 }
422};
423
433template <typename LLT, typename B_t>
434inline void llt_with_jitter(LLT& llt_B, B_t& B, double min_jitter = 1e-10,
435 double max_jitter = 1e-5) {
436 llt_B.compute(B);
437 if (llt_B.info() != Eigen::Success) {
438 double jitter_try = min_jitter;
439 for (; jitter_try < max_jitter; jitter_try *= 10) {
440 B.diagonal().array() += jitter_try;
441 llt_B.compute(B);
442 if (llt_B.info() == Eigen::Success) {
443 break;
444 }
445 }
446 if (llt_B.info() != Eigen::Success) {
447 throw std::domain_error(
448 "laplace_marginal_density: Cholesky failed after adding jitter up to "
449 + std::to_string(jitter_try));
450 }
451 }
452}
453
469 Eigen::VectorXd W_r_diag;
470
472 Eigen::VectorXd W_diag;
473
475 Eigen::LLT<Eigen::MatrixXd> llt_B;
476
477 template <typename NewtonStateT, typename CovarMat>
478 CholeskyWSolverDiag(const NewtonStateT& state, const CovarMat& covariance)
479 : W_r_diag(Eigen::VectorXd::Zero(state.b.size())), W_diag(0), llt_B() {}
499 template <typename NewtonStateT, typename LLFun, typename LLTupleArgs,
500 typename CovarMat>
501 void solve_step(NewtonStateT& state, const LLFun& ll_fun,
502 const LLTupleArgs& ll_args, const CovarMat& covariance,
503 int /*hessian_block_size*/, std::ostream* msgs) {
504 const Eigen::Index theta_size = state.b.size();
505
506 // 1. Compute diagonal Hessian
507 W_diag = laplace_likelihood::diagonal_hessian(ll_fun, state.prev().theta(),
508 ll_args, msgs);
509 for (Eigen::Index j = 0; j < W_diag.size(); j++) {
510 if (W_diag.coeff(j) < 0 || !std::isfinite(W_diag.coeff(j))) {
511 throw std::domain_error(
512 "laplace_marginal_density: Hessian matrix is not positive "
513 "definite");
514 } else {
515 W_r_diag.coeffRef(j) = std::sqrt(W_diag.coeff(j));
516 }
517 }
518
519 // 2. Formulate B = I + W_r * Sigma * W_r
520 state.B.noalias()
521 = Eigen::MatrixXd::Identity(theta_size, theta_size)
522 + W_r_diag.asDiagonal() * covariance * W_r_diag.asDiagonal();
523
524 // 3. Factorize B with jittering fallback
525 llt_with_jitter(llt_B, state.B);
526 // 4. Solve for curr.a
527 state.b.noalias() = (W_diag.array() * state.prev().theta().array()).matrix()
528 + state.prev().theta_grad();
529 auto L = llt_B.matrixL();
530 auto LT = llt_B.matrixU();
531 state.curr().a().noalias()
532 = state.b
533 - W_r_diag.asDiagonal()
534 * LT.solve(
535 L.solve(W_r_diag.cwiseProduct(covariance * state.b)));
536 }
537
542 double compute_log_determinant() const {
543 return 2.0 * llt_B.matrixLLT().diagonal().array().log().sum();
544 }
545
554 template <typename NewtonStateT>
555 auto build_result(NewtonStateT& state, double log_det) {
557 state.prev().obj() - 0.5 * log_det,
558 std::move(state).prev().theta(),
559 Eigen::SparseMatrix<double>(W_r_diag.asDiagonal()),
560 Eigen::MatrixXd(llt_B.matrixL()),
561 std::move(state).prev().a(),
562 std::move(state).prev().theta_grad(),
563 Eigen::PartialPivLU<Eigen::MatrixXd>{},
564 Eigen::MatrixXd(0, 0),
565 1};
566 }
567};
568
584 Eigen::SparseMatrix<double> W_r;
585
587 Eigen::SparseMatrix<double> W_block;
588
590 Eigen::LLT<Eigen::MatrixXd> llt_B;
591
592 template <typename NewtonStateT>
593 CholeskyWSolverBlock(const NewtonStateT& state, int hessian_block_size)
594 : W_r(state.b.size(), state.b.size()) {
595 const Eigen::Index theta_size = state.b.size();
596 W_r.reserve(Eigen::VectorXi::Constant(theta_size, hessian_block_size));
597 const Eigen::Index n_block = theta_size / hessian_block_size;
598 for (Eigen::Index ii = 0; ii < n_block; ii++) {
599 for (Eigen::Index k = 0; k < hessian_block_size; k++) {
600 for (Eigen::Index j = 0; j < hessian_block_size; j++) {
601 W_r.insert(ii * hessian_block_size + j, ii * hessian_block_size + k)
602 = 1.0;
603 }
604 }
605 }
606 W_r.makeCompressed();
607 }
608
630 template <typename NewtonStateT, typename LLFun, typename LLTupleArgs,
631 typename CovarMat>
632 void solve_step(NewtonStateT& state, const LLFun& ll_fun,
633 const LLTupleArgs& ll_args, const CovarMat& covariance,
634 int hessian_block_size, std::ostream* msgs) {
635 const Eigen::Index theta_size = state.b.size();
636 // 1. Compute block Hessian
638 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
639
640 for (Eigen::Index j = 0; j < W_block.rows(); j++) {
641 if (W_block.coeff(j, j) < 0 || !std::isfinite(W_block.coeff(j, j))) {
642 throw std::domain_error(
643 "laplace_marginal_density: Hessian matrix is not positive "
644 "definite");
645 }
646 }
647
648 // 2. Compute W_r = sqrt(W)
649 block_matrix_sqrt(W_r, W_block, hessian_block_size);
650
651 // 3. Formulate B = I + W_r * Sigma * W_r
652 state.B.noalias() = Eigen::MatrixXd::Identity(theta_size, theta_size)
653 + W_r * (covariance * W_r);
654
655 // 4. Factorize B with jittering fallback
656 llt_with_jitter(llt_B, state.B);
657
658 // 5. Solve for curr.a
659 state.b.noalias()
660 = W_block * state.prev().theta() + state.prev().theta_grad();
661 auto L = llt_B.matrixL();
662 auto LT = llt_B.matrixU();
663 state.curr().a().noalias()
664 = state.b - W_r * LT.solve(L.solve(W_r * (covariance * state.b)));
665 }
666
671 double compute_log_determinant() const {
672 return 2.0 * llt_B.matrixLLT().diagonal().array().log().sum();
673 }
674
683 template <typename NewtonStateT>
684 auto build_result(NewtonStateT& state, double log_det) {
685 return laplace_density_estimates{state.prev().obj() - 0.5 * log_det,
686 std::move(state).prev().theta(),
687 std::move(W_r),
688 Eigen::MatrixXd(llt_B.matrixL()),
689 std::move(state).prev().a(),
690 std::move(state).prev().theta_grad(),
691 Eigen::PartialPivLU<Eigen::MatrixXd>{},
692 Eigen::MatrixXd(0, 0),
693 1};
694 }
695};
696
711 Eigen::MatrixXd K_root;
712
714 Eigen::SparseMatrix<double> W_full;
715
717 Eigen::LLT<Eigen::MatrixXd> llt_B;
718
719 template <typename NewtonStateT, typename CovarMat>
720 CholeskyKSolver(const NewtonStateT& state, const CovarMat& covariance)
721 : K_root(0, 0), W_full(0, 0), llt_B() {
722 auto K_root_llt = covariance.template selfadjointView<Eigen::Lower>().llt();
723 if (K_root_llt.info() != Eigen::Success) {
724 throw std::domain_error(
725 "laplace_marginal_density: Cholesky of covariance failed at start");
726 }
727 K_root = std::move(K_root_llt.matrixL());
728 }
729
750 template <typename NewtonStateT, typename LLFun, typename LLTupleArgs,
751 typename CovarMat>
752 void solve_step(NewtonStateT& state, const LLFun& ll_fun,
753 const LLTupleArgs& ll_args, const CovarMat& covariance,
754 int hessian_block_size, std::ostream* msgs) {
755 const Eigen::Index theta_size = state.b.size();
756
757 // 1. Compute Hessian
759 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
760
761 // 2. Formulate B = I + K^T * W * K
762 state.B.noalias() = Eigen::MatrixXd::Identity(theta_size, theta_size)
763 + K_root.transpose() * (W_full * K_root);
764
765 // 3. Factorize B with jittering fallback
766 llt_with_jitter(llt_B, state.B);
767
768 // 4. Solve for curr.a
769 state.b.noalias()
770 = W_full * state.prev().theta() + state.prev().theta_grad();
771 auto L = llt_B.matrixL();
772 auto LT = llt_B.matrixU();
773 state.curr().a().noalias()
774 = K_root.transpose().template triangularView<Eigen::Upper>().solve(
775 LT.solve(L.solve(K_root.transpose() * state.b)));
776 }
777
782 double compute_log_determinant() const {
783 return 2.0 * llt_B.matrixLLT().diagonal().array().log().sum();
784 }
785
794 template <typename NewtonStateT>
795 auto build_result(NewtonStateT& state, double log_det) {
796 return laplace_density_estimates{state.prev().obj() - 0.5 * log_det,
797 std::move(state.prev().theta()),
798 std::move(W_full),
799 Eigen::MatrixXd(llt_B.matrixL()),
800 std::move(state.prev().a()),
801 std::move(state.prev().theta_grad()),
802 Eigen::PartialPivLU<Eigen::MatrixXd>{},
803 std::move(K_root),
804 2};
805 }
806};
807
821struct LUSolver {
823 Eigen::PartialPivLU<Eigen::MatrixXd> lu;
824
826 Eigen::SparseMatrix<double> W_full;
827
845 template <typename NewtonStateT, typename LLFun, typename LLTupleArgs,
846 typename CovarMat>
847 void solve_step(NewtonStateT& state, const LLFun& ll_fun,
848 const LLTupleArgs& ll_args, const CovarMat& covariance,
849 int hessian_block_size, std::ostream* msgs) {
850 const Eigen::Index theta_size = state.b.size();
851
852 // 1. Compute Hessian
854 ll_fun, state.prev().theta(), hessian_block_size, ll_args, msgs);
855
856 // 2. Factorize B = I + Sigma * W
857 lu.compute(Eigen::MatrixXd::Identity(theta_size, theta_size)
858 + covariance * W_full);
859
860 // 3. Solve for curr.a
861 state.b.noalias()
862 = W_full * state.prev().theta() + state.prev().theta_grad();
863 state.curr().a().noalias()
864 = state.b - W_full * lu.solve(covariance * state.b);
865 }
866
876 double compute_log_determinant() const {
877 return lu.matrixLU().diagonal().array().log().sum();
878 }
879
888 template <typename NewtonStateT>
889 auto build_result(NewtonStateT& state, double log_det) {
890 return laplace_density_estimates{state.prev().obj() - 0.5 * log_det,
891 std::move(state).prev().theta(),
892 std::move(W_full),
893 Eigen::MatrixXd(0, 0),
894 std::move(state).prev().a(),
895 std::move(state).prev().theta_grad(),
896 std::move(lu),
897 Eigen::MatrixXd(0, 0),
898 3};
899 }
900};
901
924template <typename SolverPolicy, typename NewtonStateT, typename OptionsT,
925 typename LLFunT, typename LLTupleArgsT, typename CovarMatT,
926 typename UpdateFun>
927inline auto run_newton_loop(SolverPolicy& solver, NewtonStateT& state,
928 const OptionsT& options, Eigen::Index& step_iter,
929 const LLFunT& ll_fun, const LLTupleArgsT& ll_args,
930 const CovarMatT& covariance, UpdateFun&& update_fun,
931 std::ostream* msgs) {
932 bool finish_update = false;
933 for (; step_iter <= options.max_num_steps; step_iter++) {
934 solver.solve_step(state, ll_fun, ll_args, covariance,
935 options.hessian_block_size, msgs);
936 if (!state.final_loop) {
937 state.wolfe_info.p_ = state.curr().a() - state.prev().a();
938 state.prev_g.noalias() = -covariance * state.prev().a()
939 + covariance * state.prev().theta_grad();
940 state.wolfe_info.init_dir_ = state.prev_g.dot(state.wolfe_info.p_);
941 // Flip direction if not ascending
942 state.wolfe_info.flip_direction();
943 auto&& scratch = state.wolfe_info.scratch_;
944 scratch.alpha() = 1.0;
945 update_fun(scratch, state.curr(), state.prev(), scratch.eval_,
946 state.wolfe_info.p_);
947 if (scratch.alpha() <= options.line_search.min_alpha) {
948 state.wolfe_status.accept_ = false;
949 finish_update = true;
950 } else if (options.line_search.max_iterations == 0) {
951 state.curr().update(scratch);
952 state.wolfe_status.accept_ = true;
953 } else {
954 Eigen::VectorXd s = scratch.a() - state.prev().a();
955 auto full_step_grad
956 = (-covariance * scratch.a() + covariance * scratch.theta_grad())
957 .eval();
958 state.curr().alpha() = barzilai_borwein_step_size(
959 s, full_step_grad, state.prev_g, state.prev().alpha(),
960 state.wolfe_status.num_backtracks_, options.line_search.min_alpha,
961 options.line_search.max_alpha);
962 state.wolfe_status = internal::wolfe_line_search(
963 state.wolfe_info, update_fun, options.line_search, msgs);
964 }
969 bool objective_converged
970 = std::abs(state.curr().obj() - state.prev().obj())
971 < options.tolerance;
972 bool search_failed = (!state.wolfe_status.accept_
973 && state.curr().obj() <= state.prev().obj());
974 finish_update = objective_converged || search_failed;
975 }
976 if (finish_update) {
977 if (!state.final_loop && state.wolfe_status.accept_) {
978 // Do one final loop with exact wolfe conditions
979 state.final_loop = true;
980 // NOTE: Swapping here so we need to swap prev and curr later
981 state.update_next_step(options);
982 continue;
983 }
984 return solver.build_result(state, solver.compute_log_determinant());
985 } else {
986 state.update_next_step(options);
987 }
988 }
989 if (msgs) {
990 (*msgs)
991 << std::string(
992 "WARNING(laplace_marginal_density): max number of iterations: ")
993 + std::to_string(options.max_num_steps) + " exceeded.";
994 }
995 return solver.build_result(state, solver.compute_log_determinant());
996}
997
1008inline void log_solver_fallback(const bool allow_fallthrough,
1009 std::ostream* msgs, std::string_view context,
1010 Eigen::Index iter,
1011 std::string_view failed_solver,
1012 std::string_view next_solver,
1013 const std::exception& e) {
1014 // Build once so we don't interleave with other logs.
1015 std::ostringstream os;
1016 std::string msg_type = allow_fallthrough ? "WARNING" : "ERROR";
1017 os << "[" << context << "] " << msg_type << ": solver fallback\n"
1018 << " " << std::left << std::setw(12) << "iteration:" << iter << "\n"
1019 << " " << std::left << std::setw(12) << "failed:" << failed_solver << "\n"
1020 << " " << std::left << std::setw(12) << "reason:" << e.what() << "\n"
1021 << " " << std::left << std::setw(12) << "action:"
1022 << "trying " << next_solver << "\n"
1023 << "note: this warning message will only be displayed once."
1024 << "\n";
1025 if (allow_fallthrough && msgs) {
1026 (*msgs) << os.str();
1027 } else {
1028 throw std::domain_error(std::string("[") + std::string(context) + "]");
1029 }
1030}
1031
1032template <bool InitTheta, typename Opts>
1033inline decltype(auto) theta_init_impl(Eigen::Index theta_size, Opts&& options) {
1034 if constexpr (InitTheta) {
1035 // If requested, use the prior mean as the initial value
1036 return std::decay_t<decltype(options)>(options).theta_0;
1037 } else {
1038 return Eigen::MatrixXd::Zero(theta_size, 1);
1039 }
1040}
1041
1060template <typename ObjFun, typename ThetaGradFun, typename Covariance,
1061 typename Options>
1062inline auto create_update_fun(ObjFun&& obj_fun, ThetaGradFun&& theta_grad_f,
1063 Covariance&& covariance, Options&& options) {
1064 auto update_step = [&covariance, &obj_fun, &theta_grad_f](
1065 auto& proposal, auto&& /* curr */, auto&& prev,
1066 auto& eval_in, auto&& p) {
1067 try {
1068 proposal.a() = prev.a() + eval_in.alpha() * p;
1069 proposal.theta().noalias() = covariance * proposal.a();
1070 proposal.theta_grad() = theta_grad_f(proposal.theta());
1071 eval_in.obj() = obj_fun(proposal.a(), proposal.theta());
1072 eval_in.dir()
1073 = (-covariance * proposal.a() + covariance * proposal.theta_grad())
1074 .dot(p);
1075 return std::isfinite(eval_in.obj()) && std::isfinite(eval_in.dir());
1076 } catch (const std::exception&) {
1077 return false;
1078 }
1079 };
1080 auto backoff = [&options](auto& eval) {
1081 eval.alpha() *= options.line_search.tau;
1082 return eval.alpha() > options.line_search.min_alpha;
1083 };
1084 return
1085 [update_step_ = std::move(update_step), backoff_ = std::move(backoff)](
1086 auto& proposal, auto&& curr, auto&& prev, auto& eval_in, auto&& p) {
1087 return internal::retry_evaluate(update_step_, proposal, curr, prev,
1088 eval_in, p, backoff_);
1089 };
1090}
1091
1144template <typename LLFun, typename LLTupleArgs, typename CovarMat,
1145 bool InitTheta,
1148 LLFun&& ll_fun, LLTupleArgs&& ll_args, CovarMat&& covariance,
1149 const laplace_options<InitTheta>& options, std::ostream* msgs) {
1150 internal::validate_laplace_options("laplace_marginal_density", options,
1151 covariance);
1152 const Eigen::Index theta_size = covariance.rows();
1153 // Wolfe optimizes over the latent 'a' space
1154 auto obj_fun = [&ll_fun, &ll_args, &msgs](const Eigen::VectorXd& a_val,
1155 auto&& theta_val) -> double {
1156 return -0.5 * a_val.dot(theta_val)
1157 + laplace_likelihood::log_likelihood(ll_fun, theta_val, ll_args,
1158 msgs);
1159 };
1160 auto theta_grad_f = [&ll_fun, &ll_args, &msgs](auto&& theta_val) {
1161 return laplace_likelihood::theta_grad(ll_fun, theta_val, ll_args, msgs);
1162 };
1163 decltype(auto) theta_init = theta_init_impl<InitTheta>(theta_size, options);
1164 internal::NewtonState state(theta_size, obj_fun, theta_grad_f, theta_init);
1165 // Start with safe step size
1166 auto update_fun = create_update_fun(
1167 std::move(obj_fun), std::move(theta_grad_f), covariance, options);
1168 Eigen::Index step_iter = 0;
1169 try {
1170 if (options.solver == 1) {
1171 if (options.hessian_block_size == 1) {
1172 CholeskyWSolverDiag solver(state, covariance);
1173 return run_newton_loop(solver, state, options, step_iter, ll_fun,
1174 ll_args, covariance, update_fun, msgs);
1175 } else {
1176 CholeskyWSolverBlock solver(state, options.hessian_block_size);
1177 return run_newton_loop(solver, state, options, step_iter, ll_fun,
1178 ll_args, covariance, update_fun, msgs);
1179 }
1180 }
1181 } catch (const std::exception& e) {
1182 const std::string solver_type
1183 = (options.hessian_block_size == 1) ? "Diagonal" : "Block";
1184 std::string failed = "solver 1 (" + solver_type + " Hessian-root Cholesky)";
1185 std::call_once(
1187 [](auto&&... args) {
1188 log_solver_fallback(std::forward<decltype(args)>(args)...);
1189 },
1190 options.allow_fallthrough, msgs, "laplace_marginal_density", step_iter,
1191 std::move(failed), "solver 2 (Covariance-root Cholesky)", e);
1192 }
1193 try {
1194 if (options.solver == 2 || options.allow_fallthrough) {
1195 CholeskyKSolver solver(state, covariance);
1196 return run_newton_loop(solver, state, options, step_iter, ll_fun, ll_args,
1197 covariance, update_fun, msgs);
1198 }
1199 } catch (const std::exception& e) {
1200 std::call_once(
1202 [](auto&&... args) {
1203 log_solver_fallback(std::forward<decltype(args)>(args)...);
1204 },
1205 options.allow_fallthrough, msgs, "laplace_marginal_density", step_iter,
1206 "solver 2 (Covariance-root Cholesky)", "solver 3 (General LU solver)",
1207 e);
1208 }
1209 if (options.solver == 3 || options.allow_fallthrough) {
1210 LUSolver solver;
1211 return run_newton_loop(solver, state, options, step_iter, ll_fun, ll_args,
1212 covariance, update_fun, msgs);
1213 }
1214 throw std::domain_error(
1215 std::string("You chose a solver (") + std::to_string(options.solver)
1216 + ") that is not valid. Please choose either 1, 2, or 3.");
1217}
1218} // namespace internal
1219} // namespace math
1220} // namespace stan
1221#endif
#define STAN_THREADS_DEF
#define unlikely(x)
int64_t size(const T &m)
Returns the size (number of the elements) of a matrix_cl or var_value<matrix_cl<T>>.
Definition size.hpp:19
(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.
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 auto tuple_to_laplace_options(Options &&ops)
void validate_laplace_options(const char *frame_name, const laplace_options< InitTheta > &options, const CovarMat &covariance)
Validates the options for the Laplace approximation.
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.
Definition constants.hpp:20
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
auto generate_laplace_options(int theta_0_size)
User function for generating laplace options tuple.
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)
Definition dot.hpp:11
constexpr bool is_inner_tuple_type_v
Checks if the N-th element of a tuple is of the same type as CheckType.
Definition is_tuple.hpp:82
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 ...
STL namespace.
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.
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)
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.
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(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(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.