Accelerators#
- group grp_Accelerators
Computing accelerated directions, such as L-BFGS.
-
template<Config Conf = DefaultConfig>
class AndersonAccel - #include <alpaqa/accelerators/anderson.hpp>
Anderson Acceleration.
Algorithm for accelerating fixed-point iterations for finding fixed points of a function \( g \), i.e. \( g(x^\star) = x^\star \), or equivalently, roots of the residual \( r(x) \triangleq g(x) - x \).
Public Types
-
using Params = AndersonAccelParams<config_t>
Public Functions
-
AndersonAccel() = default
-
inline AndersonAccel(Params params)
- Parameters:
params – Parameters.
-
inline AndersonAccel(Params params, length_t n)
- Parameters:
params – Parameters.
n – Problem dimension (size of the vectors).
-
inline void resize(length_t n)
Change the problem dimension.
Flushes the history.
- Parameters:
n – Problem dimension (size of the vectors).
-
inline void initialize(crvec g_0, crvec r_0)
Call this function on the first iteration to initialize the accelerator.
-
inline void compute(crvec gₖ, crvec rₖ, rvec xₖ_aa)
Compute the accelerated iterate \( x^k_\text{AA} \), given the function value at the current iterate \( g^k = g(x^k) \) and the corresponding residual \( r^k = g^k - x^k \).
-
inline void compute(crvec gₖ, vec &&rₖ, rvec xₖ_aa)
Compute the accelerated iterate \( x^k_\text{AA} \), given the function value at the current iterate \( g^k = g(x^k) \) and the corresponding residual \( r^k = g^k - x^k \).
-
inline void reset()
Reset the accelerator (but keep the last function value and residual, so calling initialize is not necessary).
-
inline length_t n() const
Get the problem dimension.
-
inline length_t history() const
Get the maximum number of stored columns.
-
inline length_t current_history() const
Get the number of columns currently stored in the buffer.
-
inline const Params &get_params() const
Get the parameters.
-
inline std::string get_name() const
-
inline void scale_R(real_t scal)
Scale the factorization.
-
inline const LimitedMemoryQR<config_t> &get_QR() const
For testing purposes.
-
using Params = AndersonAccelParams<config_t>
-
template<Config Conf = DefaultConfig>
class LBFGS - #include <alpaqa/accelerators/lbfgs.hpp>
Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm.
Public Types
-
enum class Sign
The sign of the vectors \( p \) passed to the update method.
Values:
-
enumerator Positive
\( p \sim \nabla \psi(x) \)
-
enumerator Negative
\( p \sim -\nabla \psi(x) \)
-
enumerator Positive
-
using Params = LBFGSParams<config_t>
Public Functions
-
LBFGS() = default
-
inline LBFGS(Params params)
-
bool update_sy(crvec s, crvec y, real_t pₙₑₓₜᵀpₙₑₓₜ, bool forced = false)
Update the inverse Hessian approximation using the new vectors sₖ = xₙₑₓₜ - xₖ and yₖ = pₙₑₓₜ - pₖ.
-
bool update_sy_impl(const auto &s, const auto &y, real_t pₙₑₓₜᵀpₙₑₓₜ, bool forced = false)
See also
-
bool update(crvec xₖ, crvec xₙₑₓₜ, crvec pₖ, crvec pₙₑₓₜ, Sign sign = Sign::Positive, bool forced = false)
Update the inverse Hessian approximation using the new vectors xₙₑₓₜ and pₙₑₓₜ.
-
bool apply(rvec q, real_t γ = -1) const
Apply the inverse Hessian approximation to the given vector q.
Initial inverse Hessian approximation is set to \( H_0 = \gamma I \). If
γ
is negative, \( H_0 = \frac{s^\top y}{y^\top y} I \).
-
bool apply_masked(rvec q, real_t γ, crindexvec J) const
Apply the inverse Hessian approximation to the given vector q, applying only the columns and rows of the Hessian in the index set J.
-
bool apply_masked(rvec q, real_t γ, const std::vector<index_t> &J) const
Apply the inverse Hessian approximation to the given vector q, applying only the columns and rows of the Hessian in the index set J.
-
bool apply_masked_impl(rvec q, real_t γ, const auto &J) const
Apply the inverse Hessian approximation to the given vector q, applying only the columns and rows of the Hessian in the index set J.
-
void reset()
Throw away the approximation and all previous vectors s and y.
-
void scale_y(real_t factor)
Scale the stored y vectors by the given factor.
-
inline std::string get_name() const
Get a string identifier for this accelerator.
-
inline const Params &get_params() const
Get the parameters.
-
inline length_t n() const
Get the size of the s and y vectors in the buffer.
-
inline length_t history() const
Get the number of previous vectors s and y stored in the buffer.
-
inline index_t succ(index_t i) const
Get the next index in the circular buffer of previous s and y vectors.
-
inline index_t pred(index_t i) const
Get the previous index in the circular buffer of s and y vectors.
-
inline length_t current_history() const
Get the number of previous s and y vectors currently stored in the buffer.
-
inline auto s(index_t i)
-
inline auto s(index_t i) const
-
inline auto y(index_t i)
-
inline auto y(index_t i) const
-
template<class F>
inline void foreach_fwd(const F &fun) const Iterate over the indices in the history buffer, oldest first.
-
template<class F>
inline void foreach_rev(const F &fun) const Iterate over the indices in the history buffer, newest first.
-
enum class Sign
-
template<Config Conf = DefaultConfig>