quala 0.0.1a1
Quasi-Newton and other accelerators
Namespaces | Classes | Typedefs | Functions | Variables
quala Namespace Reference

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ₖ, 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...
 

Class Documentation

◆ quala::AndersonAccelParams

struct quala::AndersonAccelParams
+ Collaboration diagram for 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).

◆ quala::BroydenGoodParams

struct quala::BroydenGoodParams
+ Collaboration diagram for 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.

◆ quala::LBFGSParams

struct quala::LBFGSParams
+ Collaboration diagram for 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 \]

See also
https://epubs.siam.org/doi/10.1137/S1052623499354242
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 \).

Typedef Documentation

◆ real_t

using real_t = double

Default floating point type.

Definition at line 8 of file vec.hpp.

◆ realvec

using realvec = Eigen::Matrix<real_t, Eigen::Dynamic, 1>

Default type for floating point vectors.

Definition at line 10 of file vec.hpp.

◆ realmat

using realmat = Eigen::Matrix<real_t, Eigen::Dynamic, Eigen::Dynamic>

Default type for floating point matrices.

Definition at line 12 of file vec.hpp.

◆ vec

using vec = realvec

Default type for vectors.

Definition at line 14 of file vec.hpp.

◆ rvec

using rvec = Eigen::Ref<vec>

Default type for mutable references to vectors.

Definition at line 16 of file vec.hpp.

◆ crvec

using crvec = Eigen::Ref<const vec>

Default type for immutable references to vectors.

Definition at line 18 of file vec.hpp.

◆ mat

using mat = realmat

Default type for matrices.

Definition at line 20 of file vec.hpp.

◆ rmat

using rmat = Eigen::Ref<mat>

Default type for mutable references to matrices.

Definition at line 22 of file vec.hpp.

◆ crmat

using crmat = Eigen::Ref<const mat>

Default type for immutable references to matrices.

Definition at line 24 of file vec.hpp.

◆ index_t

using index_t = Eigen::Index

Default type for vector indices.

Definition at line 27 of file vec.hpp.

◆ length_t

using length_t = index_t

Default type for vector sizes.

Definition at line 29 of file vec.hpp.

◆ idvec

using idvec = Eigen::Matrix<index_t, Eigen::Dynamic, 1>

Type for a vector of indices.

Definition at line 32 of file vec.hpp.

◆ ridvec

using ridvec = Eigen::Ref<idvec>

Mutable reference to vector indices.

Definition at line 34 of file vec.hpp.

◆ cridvec

using cridvec = Eigen::Ref<const idvec>

Immutable reference to vector indices.

Definition at line 36 of file vec.hpp.

◆ anymat

using anymat = Eigen::MatrixBase<Derived>

Generic type for vector and matrix arguments.

Definition at line 40 of file vec.hpp.

Function Documentation

◆ minimize_update_anderson()

void minimize_update_anderson ( LimitedMemoryQR qr,
rmat  G,
crvec  rₖ,
crvec  rₗₐₛₜ,
crvec  gₖ,
rvec  γ_LS,
rvec  xₖ_aa 
)
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} \mathcal{R}_k &= \begin{pmatrix} \Delta r_{k-m_k} & \dots & \Delta r_{k-1} \end{pmatrix} \in \mathbb{R}^{n\times m_k} \\ \Delta r_i &= r_{i+1} - r_i \\ r_i &= g_i - x_i \\ g_i &= g(x_i) \\ \DeclareMathOperator*{\argmin}{argmin} \gamma_\text{LS} &= \argmin_\gamma \left\| \mathcal{R}_k \gamma - r_k \right\|_2 \\ \alpha_i &= \begin{cases} \gamma_\text{LS}[0] & i = 0 \\ \gamma_\text{LS}[i] - \gamma_\text{LS}[i-1] & 0 < i < m_k \\ 1 - \gamma_\text{LS}[m_k - 1] & i = m_k \end{cases} \\ x_{k+1} &= \sum_{i=0}^{m_k} \alpha_i\,g_i \end{aligned} \]

Definition at line 30 of file anderson-helpers.hpp.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Variable Documentation

◆ inf

constexpr real_t inf = std::numeric_limits<real_t>::infinity()
constexpr

\( \infty \)

Definition at line 43 of file vec.hpp.

◆ NaN

constexpr real_t NaN = std::numeric_limits<real_t>::quiet_NaN()
constexpr

Not a number.

Definition at line 45 of file vec.hpp.