quala 0.0.1a1
Quasi-Newton and other accelerators
anderson-acceleration.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <stdexcept>
5
6namespace quala {
7
8/// Parameters for the @ref AndersonAccel class.
10 /// Length of the history to keep (the number of columns in the QR
11 /// factorization).
12 /// If this number is greater than the problem dimension, the memory is set
13 /// to the problem dimension (otherwise the system is underdetermined).
15};
16
17/**
18 * Anderson Acceleration.
19 *
20 * Algorithm for accelerating fixed-point iterations for finding fixed points
21 * of a function @f$ g @f$, i.e. @f$ g(x^\star) = x^\star @f$, or equivalently,
22 * roots of the residual @f$ r(x) \triangleq g(x) - x @f$.
23 *
24 * @todo Condition estimation of the QR factorization.
25 *
26 * @ingroup accelerators-grp
27 */
29 public:
31
32 /// @param params
33 /// Parameters.
35 /// @param params
36 /// Parameters.
37 /// @param n
38 /// Problem dimension (size of the vectors).
40
41 /// Change the problem dimension. Flushes the history.
42 /// @param n
43 /// Problem dimension (size of the vectors).
45 length_t m_AA = std::min(n, params.memory); // TODO: support m > n?
46 qr.resize(n, m_AA);
47 G.resize(n, m_AA);
48 rₗₐₛₜ.resize(n);
49 γ_LS.resize(m_AA);
50 initialized = false;
51 }
52
53 /// Call this function on the first iteration to initialize the accelerator.
54 void initialize(crvec g_0, vec r_0) {
55 assert(g_0.size() == vec::Index(n()));
56 assert(r_0.size() == vec::Index(n()));
57 G.col(0) = g_0;
58 rₗₐₛₜ = std::move(r_0);
59 qr.reset();
60 initialized = true;
61 }
62
63 /// Compute the accelerated iterate @f$ x^k_\text{AA} @f$, given the
64 /// function value at the current iterate @f$ g^k = g(x^k) @f$ and the
65 /// corresponding residual @f$ r^k = g^k - x^k @f$.
66 void compute(crvec gₖ, crvec rₖ, rvec xₖ_aa) {
67 if (!initialized)
68 throw std::logic_error("AndersonAccel::compute() called before "
69 "AndersonAccel::initialize()");
71 rₖ, rₗₐₛₜ, gₖ, // in
72 γ_LS, xₖ_aa); // out
73 rₗₐₛₜ = rₖ;
74 }
75 /// @copydoc compute(crvec, crvec, rvec)
76 void compute(crvec gₖ, vec &&rₖ, rvec xₖ_aa) {
77 if (!initialized)
78 throw std::logic_error("AndersonAccel::compute() called before "
79 "AndersonAccel::initialize()");
81 rₖ, rₗₐₛₜ, gₖ, // in
82 γ_LS, xₖ_aa); // out
83 rₗₐₛₜ = std::move(rₖ);
84 }
85
86 /// Reset the accelerator (but keep the last function value and residual, so
87 /// calling @ref initialize is not necessary).
88 void reset() {
89 index_t newest_g_idx = qr.ring_tail();
90 if (newest_g_idx != 0)
91 G.col(0) = G.col(newest_g_idx);
92 qr.reset();
93 }
94
95 /// Get the problem dimension.
96 length_t n() const { return qr.n(); }
97 /// Get the maximum number of stored columns.
98 length_t history() const { return qr.m(); }
99 /// Get the number of columns currently stored in the buffer.
101
102 /// Get the parameters.
103 const Params &get_params() const { return params; }
104
105 private:
111 bool initialized = false;
112};
113
114} // namespace quala
Anderson Acceleration.
AndersonAccel(Params params, length_t n)
length_t current_history() const
Get the number of columns currently stored in the buffer.
length_t n() const
Get the problem dimension.
void compute(crvec gₖ, crvec rₖ, rvec xₖ_aa)
Compute the accelerated iterate , given the function value at the current iterate and the correspond...
void initialize(crvec g_0, vec r_0)
Call this function on the first iteration to initialize the accelerator.
const Params & get_params() const
Get the parameters.
length_t history() const
Get the maximum number of stored columns.
void resize(length_t n)
Change the problem dimension.
void reset()
Reset the accelerator (but keep the last function value and residual, so calling initialize is not ne...
void compute(crvec gₖ, vec &&rₖ, rvec xₖ_aa)
Compute the accelerated iterate , given the function value at the current iterate and the correspond...
Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.
length_t current_history() const
Get the number of columns currently stored in the buffer.
index_t ring_tail() const
Get the tail index of the circular buffer (points to one past the most recent element).
void resize(length_t n, length_t m)
Re-allocate storage for a problem with a different size.
void reset()
Reset all indices, clearing the Q and R matrices.
Eigen::Ref< const vec > crvec
Default type for immutable references to vectors.
Definition: vec.hpp:18
length_t memory
Length of the history to keep (the number of columns in the QR factorization).
realmat mat
Default type for matrices.
Definition: vec.hpp:20
index_t length_t
Default type for vector sizes.
Definition: vec.hpp:29
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):
realvec vec
Default type for vectors.
Definition: vec.hpp:14
Eigen::Index index_t
Default type for vector indices.
Definition: vec.hpp:27
Eigen::Ref< vec > rvec
Default type for mutable references to vectors.
Definition: vec.hpp:16
Parameters for the AndersonAccel class.