Namespaces | |
namespace | vec_util |
Classes | |
class | AndersonAccel |
Anderson Acceleration. More... | |
struct | AndersonAccelParams |
Parameters for the AndersonAccel class. More... | |
class | BroydenGood |
Broyden's “Good” method for solving systems of nonlinear equations \( F(x) = 0 \). More... | |
struct | BroydenGoodParams |
Parameters for the BroydenGood class. More... | |
struct | BroydenStorage |
Layout: More... | |
struct | CircularIndexIterator |
struct | CircularIndices |
class | CircularRange |
class | LBFGS |
Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm. More... | |
struct | LBFGSParams |
Parameters for the LBFGS class. More... | |
struct | LBFGSStorage |
Layout: More... | |
class | LimitedMemoryQR |
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization. More... | |
class | MaxHistory |
struct | ReverseCircularIndexIterator |
class | ReverseCircularRange |
class | vec_allocator |
Typedefs | |
using | real_t = double |
Default floating point type. More... | |
using | realvec = Eigen::Matrix< real_t, Eigen::Dynamic, 1 > |
Default type for floating point vectors. More... | |
using | realmat = Eigen::Matrix< real_t, Eigen::Dynamic, Eigen::Dynamic > |
Default type for floating point matrices. More... | |
using | vec = realvec |
Default type for vectors. More... | |
using | rvec = Eigen::Ref< vec > |
Default type for mutable references to vectors. More... | |
using | crvec = Eigen::Ref< const vec > |
Default type for immutable references to vectors. More... | |
using | mat = realmat |
Default type for matrices. More... | |
using | rmat = Eigen::Ref< mat > |
Default type for mutable references to matrices. More... | |
using | crmat = Eigen::Ref< const mat > |
Default type for immutable references to matrices. More... | |
using | index_t = Eigen::Index |
Default type for vector indices. More... | |
using | length_t = index_t |
Default type for vector sizes. More... | |
using | idvec = Eigen::Matrix< index_t, Eigen::Dynamic, 1 > |
Type for a vector of indices. More... | |
using | ridvec = Eigen::Ref< idvec > |
Mutable reference to vector indices. More... | |
using | cridvec = Eigen::Ref< const idvec > |
Immutable reference to vector indices. More... | |
template<class Derived > | |
using | anymat = Eigen::MatrixBase< Derived > |
Generic type for vector and matrix arguments. More... | |
Functions | |
void | minimize_update_anderson (LimitedMemoryQR &qr, rmat G̃, crvec rₖ, crvec rₗₐₛₜ, crvec gₖ, real_t min_div, rvec γ_LS, rvec xₖ_aa) |
Solve one step of Anderson acceleration to find a fixed point of a function g(x): More... | |
Variables | |
constexpr real_t | inf = std::numeric_limits<real_t>::infinity() |
\( \infty \) More... | |
constexpr real_t | NaN = std::numeric_limits<real_t>::quiet_NaN() |
Not a number. More... | |
struct quala::AndersonAccelParams |
Class Members | ||
---|---|---|
length_t | memory = 10 |
Length of the history to keep (the number of columns in the QR factorization). If this number is greater than the problem dimension, the memory is set to the problem dimension (otherwise the system is underdetermined). |
real_t | min_div = 1e2 * std::numeric_limits<real_t>::epsilon() | Minimum divisor when solving close to singular systems, scaled by the maximum eigenvalue of R. |
struct quala::BroydenGoodParams |
Class Members | ||
---|---|---|
length_t | memory = 10 | Length of the history to keep. |
real_t | min_div_abs = 1e-32 | Reject update if \( s^\top Hy \le \text{min_div_fac} \). |
bool | force_pos_def = false | If set to true, the inverse Jacobian estimate should remain definite. |
bool | restarted = true |
If set to true, the buffer is cleared after memory iterations. If set to false, a circular buffer is used, replacing the oldest pair of vectors, as a limited-memory type algorithm. However, since each vector s̃ depends on the previous vectors, this algorithm is not theoretically correct, because the s̃ vectors are not updated when the old history is overwritten. |
real_t | powell_damping_factor = 0 | Powell's trick, damping, prevents nonsingularity. |
real_t | min_stepsize = 1e-10 |
Minimum automatic step size. If \( \frac{s^\top y}{y^\top y} \) is smaller than this setting, use this as the step size instead. |
struct quala::LBFGSParams |
Class Members | ||
---|---|---|
length_t | memory = 10 | Length of the history to keep. |
real_t | min_div_fac = 1e-10 | Reject update if \( y^\top s \le \text{min_div_fac} \cdot s^\top s \). |
real_t | min_abs_s = 1e-32 | Reject update if \( s^\top s \le \text{min_abs_s} \). |
CBFGSParams | cbfgs |
Parameters in the cautious BFGS update condition. \[ \frac{y^\top s}{s^\top s} \ge \epsilon \| g \|^\alpha \] |
bool | force_pos_def = true |
If set to true, the inverse Hessian estimate should remain definite, i.e. a check is performed that rejects the update if \( y^\top s \le \text{min_div_fac} \cdot s^\top s \). If set to false, just try to prevent a singular Hessian by rejecting the update if \( \left| y^\top s \right| \le \text{min_div_fac} \cdot s^\top s \). |
using index_t = Eigen::Index |
using anymat = Eigen::MatrixBase<Derived> |
|
inline |
Solve one step of Anderson acceleration to find a fixed point of a function g(x):
\( g(x^\star) - x^\star = 0 \)
Updates the QR factorization of \( \mathcal{R}_k = QR \), solves the least squares problem to find \( \gamma_\text{LS} \), computes the next iterate \( x_{k+1} \), and stores the current function value \( g_k \) in the matrix \( G \), which is used as a circular buffer.
\[ \begin{aligned} \def\gammaLS{\gamma_\text{LS}} m_k &= \min \{k, m\} \\ g_i &= g(x_i) \\ r_i &= r(x_i) g_i - x_i \\ \Delta r_i &= r_i - r_{i-1} \\ \mathcal{R}_k &= \begin{pmatrix} \Delta r_{k-m_k+1} & \dots & \Delta r_k \end{pmatrix} \in \R^{n\times m_k} \\ \gammaLS &= \argmin_{\gamma \in \R^{m_k}}\; \norm{\mathcal{R}_k \gamma - r_k}_2 \\ \alpha_i &= \begin{cases} \gammaLS[0] & i = 0 \\ \gammaLS[i] - \gammaLS[i-1] & 0 \lt i \lt m_k \\ 1 - \gammaLS[m_k - 1] & i = m_k \end{cases} \\ \tilde G_k &= \begin{pmatrix} g_{k - m_k} & \dots & g_{k-1} \end{pmatrix} \\ G_k &= \begin{pmatrix} g_{k - m_k} & \dots & g_{k} \end{pmatrix} \\ &= \begin{pmatrix} \tilde G_k & g_{k} \end{pmatrix} \\ x_{k+1} &= \sum_{i=0}^{m_k} \alpha_i\,g_{k - m_k + i} \\ &= G_k \alpha \\ \end{aligned} \]
[in,out] | qr | QR factorization of \( \mathcal{R}_k \) |
[in,out] | G̃ | Matrix of previous function values \( \tilde G_k \) (stored as ring buffer with the same indices as qr ) |
[in] | rₖ | Current residual \( r_k \) |
[in] | rₗₐₛₜ | Previous residual \( r_{k-1} \) |
[in] | gₖ | Current function value \( g_k \) |
[in] | min_div | Minimum divisor when solving close to singular systems, scaled by the maximum eigenvalue of R |
[out] | γ_LS | Solution to the least squares system |
[out] | xₖ_aa | Next Anderson iterate |
Definition at line 36 of file anderson-helpers.hpp.