|
| 1 | +#ifndef STAN_MATH_MIX_FUNCTOR_BARZILAI_BORWEIN_STEP_SIZE_HPP |
| 2 | +#define STAN_MATH_MIX_FUNCTOR_BARZILAI_BORWEIN_STEP_SIZE_HPP |
| 3 | +#include <stan/math/prim/fun/Eigen.hpp> |
| 4 | +#include <algorithm> |
| 5 | +#include <numeric> |
| 6 | +#include <cmath> |
| 7 | + |
| 8 | +namespace stan::math::internal { |
| 9 | +/** |
| 10 | + * @brief Curvature-aware Barzilai–Borwein (BB) step length with robust |
| 11 | + * safeguards. |
| 12 | + * |
| 13 | + * Given successive parameter displacements \f$s = x_k - x_{k-1}\f$ and |
| 14 | + * gradients \f$g_k\f$, \f$g_{k-1}\f$, this routine forms |
| 15 | + * \f$y = g_k - g_{k-1}\f$ and computes the two classical BB candidates |
| 16 | + * |
| 17 | + * \f{align*}{ |
| 18 | + * \alpha_{\text{BB1}} &= \frac{\langle s,s\rangle}{\langle s,y\rangle},\\ |
| 19 | + * \alpha_{\text{BB2}} &= \frac{\langle s,y\rangle}{\langle y,y\rangle}, |
| 20 | + * \f} |
| 21 | + * |
| 22 | + * then chooses between them using the **spectral cosine** |
| 23 | + * \f$r = \cos^2\!\angle(s,y) = \dfrac{\langle s,y\rangle^2} |
| 24 | + * {\langle s,s\rangle\,\langle |
| 25 | + * y,y\rangle}\in[0,1]\f$: |
| 26 | + * |
| 27 | + * - if \f$r > 0.9\f$ (well-aligned curvature) and the previous line search |
| 28 | + * did **≤ 1** backtrack, prefer the “long” step \f$\alpha_{\text{BB1}}\f$; |
| 29 | + * - if \f$0.1 \le r \le 0.9\f$, take the neutral geometric mean |
| 30 | + * \f$\sqrt{\alpha_{\text{BB1}}\alpha_{\text{BB2}}}\f$; |
| 31 | + * - otherwise default to the “short” step \f$\alpha_{\text{BB2}}\f$. |
| 32 | + * |
| 33 | + * All candidates are clamped into \f$[\text{min\_alpha},\,\text{max\_alpha}]\f$ |
| 34 | + * and must be finite and positive. |
| 35 | + * If the curvature scalars are ill-posed (non-finite or too small), |
| 36 | + * \f$\langle s,y\rangle \le \varepsilon\f$, or if `last_backtracks == 99` |
| 37 | + * (explicitly disabling BB for this iteration), the function falls back to a |
| 38 | + * **safe** step: |
| 39 | + * use `prev_step` when finite and positive, otherwise \f$1.0\f$, then clamp to |
| 40 | + * \f$[\text{min\_alpha},\,\text{max\_alpha}]\f$. |
| 41 | + * |
| 42 | + * @param s Displacement between consecutive iterates |
| 43 | + * (\f$s = x_k - x_{k-1}\f$). |
| 44 | + * @param g_curr Gradient at the current iterate \f$g_k\f$. |
| 45 | + * @param g_prev Gradient at the previous iterate \f$g_{k-1}\f$. |
| 46 | + * @param prev_step Previously accepted step length (used by the fallback). |
| 47 | + * @param last_backtracks |
| 48 | + * Number of backtracking contractions performed by the most |
| 49 | + * recent line search; set to 99 to force the safe fallback. |
| 50 | + * @param min_alpha Lower bound for the returned step length. |
| 51 | + * @param max_alpha Upper bound for the returned step length. |
| 52 | + * |
| 53 | + * @return A finite, positive BB-style step length \f$\alpha \in |
| 54 | + * [\text{min\_alpha},\,\text{max\_alpha}]\f$ suitable for seeding a |
| 55 | + * line search or as a spectral preconditioner scale. |
| 56 | + * |
| 57 | + * @note Uses \f$\varepsilon=10^{-16}\f$ to guard against division by very |
| 58 | + * small curvature terms, and applies `std::abs` to BB ratios to avoid |
| 59 | + * negative steps; descent is enforced by the line search. |
| 60 | + * @warning The vectors must have identical size. Non-finite inputs yield the |
| 61 | + * safe fallback. |
| 62 | + */ |
| 63 | +inline double barzilai_borwein_step_size(const Eigen::VectorXd& s, |
| 64 | + const Eigen::VectorXd& g_curr, |
| 65 | + const Eigen::VectorXd& g_prev, |
| 66 | + double prev_step, int last_backtracks, |
| 67 | + double min_alpha, double max_alpha) { |
| 68 | + // Fallbacks |
| 69 | + auto safe_fallback = [&]() -> double { |
| 70 | + double a = std::clamp( |
| 71 | + prev_step > 0.0 && std::isfinite(prev_step) ? prev_step : 1.0, |
| 72 | + min_alpha, max_alpha); |
| 73 | + return a; |
| 74 | + }; |
| 75 | + |
| 76 | + const Eigen::VectorXd y = g_curr - g_prev; |
| 77 | + const double sty = s.dot(y); |
| 78 | + const double sts = s.squaredNorm(); |
| 79 | + const double yty = y.squaredNorm(); |
| 80 | + |
| 81 | + // Basic validity checks |
| 82 | + constexpr double eps = 1e-16; |
| 83 | + if (!(std::isfinite(sty) && std::isfinite(sts) && std::isfinite(yty)) |
| 84 | + || sts <= eps || yty <= eps || sty <= eps || last_backtracks == -1) { |
| 85 | + return safe_fallback(); |
| 86 | + } |
| 87 | + |
| 88 | + // BB candidates |
| 89 | + double alpha_bb1 = std::clamp(std::abs(sts / sty), min_alpha, max_alpha); |
| 90 | + double alpha_bb2 = std::clamp(std::abs(sty / yty), min_alpha, max_alpha); |
| 91 | + |
| 92 | + // Safeguard candidates |
| 93 | + if (!std::isfinite(alpha_bb1) || !std::isfinite(alpha_bb2) || alpha_bb1 <= 0.0 |
| 94 | + || alpha_bb2 <= 0.0) { |
| 95 | + return safe_fallback(); |
| 96 | + } |
| 97 | + |
| 98 | + // Spectral cosine r = cos^2(angle(s, y)) in [0,1] |
| 99 | + const double r = (sty * sty) / (sts * yty); |
| 100 | + |
| 101 | + // Heuristic thresholds (robust defaults) |
| 102 | + constexpr double kLoose = 0.9; // "nice" curvature |
| 103 | + constexpr double kTight = 0.1; // "dodgy" curvature |
| 104 | + |
| 105 | + double alpha0 = alpha_bb2; // default to short BB for robustness |
| 106 | + if (r > kLoose && last_backtracks <= 1) { |
| 107 | + // Spectrum looks friendly and line search was not harsh -> try long BB |
| 108 | + alpha0 = alpha_bb1; |
| 109 | + } else if (r >= kTight && r <= kLoose) { |
| 110 | + // Neither clearly friendly nor clearly dodgy -> neutral middle |
| 111 | + alpha0 = std::sqrt(alpha_bb1 * alpha_bb2); |
| 112 | + } // else keep alpha_bb2 |
| 113 | + |
| 114 | + // Clip to user bounds |
| 115 | + alpha0 = std::clamp(alpha0, min_alpha, max_alpha); |
| 116 | + |
| 117 | + if (!std::isfinite(alpha0) || alpha0 <= 0.0) { |
| 118 | + return safe_fallback(); |
| 119 | + } |
| 120 | + return alpha0; |
| 121 | +} |
| 122 | + |
| 123 | +} // namespace stan::math::internal |
| 124 | +#endif |
0 commit comments