quala 0.0.1a0
Quasi-Newton and other accelerators
broyden-good.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <quala/util/vec.hpp>
4
5namespace quala {
6
7/// Layout:
8/// ~~~
9/// ┌───── 2 m + 1 ─────┐
10/// ┌ ┌───┬───┬───┬───┬───┐
11/// │ │ │ │ │ │ │
12/// n │ │ s │ s̃ │ s │ s̃ │ w │
13/// │ │ │ │ │ │ │
14/// └ └───┴───┴───┴───┴───┘
15/// ~~~
17 /// Re-allocate storage for a problem with a different size.
19
20 /// Get the size of the s and s̃ vectors in the buffer.
21 length_t n() const { return sto.rows(); }
22 /// Get the number of previous vectors s and s̃ stored in the buffer.
23 length_t history() const { return (sto.cols() - 1) / 2; }
24
25 auto s(index_t i) { return sto.col(2 * i); }
26 auto s(index_t i) const { return sto.col(2 * i); }
27 auto (index_t i) { return sto.col(2 * i + 1); }
28 auto (index_t i) const { return sto.col(2 * i + 1); }
29 auto work() { return sto.col(2 * history()); }
30 auto work() const { return sto.col(2 * history()); }
31
32 using storage_t =
33 Eigen::Matrix<real_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor>;
35};
36
38 sto.resize(n, history * 2 + 1);
39}
40
41/// Parameters for the @ref BroydenGood class.
43 /// Length of the history to keep.
45 /// Reject update if @f$ s^\top Hy \le \text{min_div_fac} @f$.
47 /// If set to true, the inverse Jacobian estimate should remain definite.
48 bool force_pos_def = false;
49 /// If set to true, the buffer is cleared after @p memory iterations. If
50 /// set to false, a circular buffer is used, replacing the oldest pair of
51 /// vectors, as a limited-memory type algorithm. However, since each vector
52 /// s̃ depends on the previous vectors, this algorithm is not theoretically
53 /// correct, because the s̃ vectors are not updated when the old history is
54 /// overwritten.
55 bool restarted = true;
56};
57
58/**
59 * Broyden's “Good” method for solving systems of nonlinear equations
60 * @f$ F(x) = 0 @f$.
61 *
62 * @f[ \begin{aligned}
63 * B_{k+1} &= B_k + \frac{\left(y_k - B_k s_k\right) s_k^\top}
64 * {s_k^\top s_k} \\
65 * H_{k+1} &= H_k + \frac{\left(s_k - H_k y_k\right)
66 * \left(s_k^\top H_k\right)}
67 * {s_k^\top H_k y_k} \\
68 * s_k &\triangleq x_{k+1} - x_k \\
69 * y_k &\triangleq F(x_{k+1}) - F(x_k) \\
70 * \end{aligned} @f]
71 * Where @f$ B_k @f$ approximates the Jacobian of @f$ F(x_k) @f$ and
72 * @f$ H_k \triangleq B_k^{-1} @f$.
73 *
74 * @todo Damping.
75 *
76 * @ingroup accelerators-grp
77 */
79 public:
81
84
85 /// Update the inverse Jacobian approximation using the new vectors
86 /// sₖ = xₖ₊₁ - xₖ and yₖ = pₖ₊₁ - pₖ.
87 template <class VecS, class VecY>
88 bool update_sy(const anymat<VecS> &s, const anymat<VecY> &y,
89 bool forced = false);
90
91 /// Update the inverse Jacobian approximation using the new vectors xₖ₊₁
92 /// and pₖ₊₁.
93 bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁,
94 bool forced = false);
95
96 /// Apply the inverse Jacobian approximation to the given vector q, i.e.
97 /// @f$ q \leftarrow H_k q @f$.
98 /// Initial inverse Jacobian approximation is set to @f$ H_0 = I @f$.
99 bool apply(rvec q);
100
101 /// Throw away the approximation and all previous vectors s and y.
102 void reset();
103 /// Re-allocate storage for a problem with a different size. Causes
104 /// a @ref reset.
105 void resize(length_t n);
106
107 /// Get the parameters.
108 const Params &get_params() const { return params; }
109
110 /// Get the size of the s and y vectors in the buffer.
111 length_t n() const { return sto.n(); }
112 /// Get the number of previous vectors s and y stored in the buffer.
113 length_t history() const { return sto.history(); }
114 /// Get the next index in the circular buffer of previous s and y vectors.
115 index_t succ(index_t i) const { return i + 1 < history() ? i + 1 : 0; }
116 /// Get the previous index in the circular buffer of s and y vectors.
117 index_t pred(index_t i) const { return i > 0 ? i - 1 : history() - 1; }
118 /// Get the number of previous s and y vectors currently stored in the
119 /// buffer.
120 length_t current_history() const { return full ? history() : idx; }
121
122 auto s(index_t i) { return sto.s(i); }
123 auto s(index_t i) const { return sto.s(i); }
124 auto (index_t i) { return sto.(i); }
125 auto (index_t i) const { return sto.(i); }
126
127 /// Iterate over the indices in the history buffer, oldest first.
128 template <class F>
129 void foreach_fwd(const F &fun) const {
130 if (full)
131 for (index_t i = idx; i < history(); ++i)
132 fun(i);
133 if (idx)
134 for (index_t i = 0; i < idx; ++i)
135 fun(i);
136 }
137
138 /// Iterate over the indices in the history buffer, newest first.
139 template <class F>
140 void foreach_rev(const F &fun) const {
141 if (idx)
142 for (index_t i = idx; i-- > 0;)
143 fun(i);
144 if (full)
145 for (index_t i = history(); i-- > idx;)
146 fun(i);
147 }
148
149 private:
152 bool full = false;
154};
155
156inline void BroydenGood::reset() {
157 idx = 0;
158 full = false;
159}
160
162 if (params.memory < 1)
163 throw std::invalid_argument("BroydenGood::Params::memory must be >= 1");
165 reset();
166}
167
168template <class VecS, class VecY>
170 bool forced) {
171 // Restart if the buffer is full
172 if (full && params.restarted) {
173 full = false;
174 assert(idx == 0);
175 }
176
177 auto &&r = sto.work();
178 // Compute r = r₍ₘ₋₁₎ = Hₖ yₖ
179 r = yₖ; // r₍₋₁₎ = yₖ
180 foreach_fwd([&](index_t i) {
181 r += (i) * r.dot(s(i)); // r₍ᵢ₎ = r₍ᵢ₋₁₎ + s̃₍ᵢ₎〈r₍ᵢ₋₁₎, s₍ᵢ₎〉
182 });
183 const real_t sᵀHy = sₖ.dot(r);
184 const real_t a_sᵀHy = params.force_pos_def ? sᵀHy : std::abs(sᵀHy);
185
186 // Check if update is accepted
187 if (!forced && a_sᵀHy < params.min_div_abs)
188 return false;
189
190 // Store the new vectors
191 sto.s(idx) = sₖ;
192 sto.(idx) = (1 / sᵀHy) * (sₖ - r);
193
194 // Increment the index in the circular buffer
195 idx = succ(idx);
196 full |= idx == 0;
197
198 return true;
199}
200
201inline bool BroydenGood::apply(rvec q) {
202 // Only apply if we have previous vectors s and y
203 if (idx == 0 && not full)
204 return false;
205
206 // Compute q = q₍ₘ₋₁₎ = Hₖ q
207 // q₍₋₁₎ = q
208 foreach_fwd([&](index_t i) {
209 q += (i) * q.dot(s(i)); // q₍ᵢ₎ = q₍ᵢ₋₁₎ + s̃₍ᵢ₎〈q₍ᵢ₋₁₎, s₍ᵢ₎〉
210 });
211
212 return true;
213}
214
215inline bool BroydenGood::update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁,
216 bool forced) {
217 const auto s = xₖ₊₁ - xₖ;
218 const auto y = pₖ₊₁ - pₖ;
219 return update_sy(s, y, forced);
220}
221
222} // namespace quala
Broyden's “Good” method for solving systems of nonlinear equations .
auto s̃(index_t i) const
void foreach_rev(const F &fun) const
Iterate over the indices in the history buffer, newest first.
bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, bool forced=false)
Update the inverse Jacobian approximation using the new vectors xₖ₊₁ and pₖ₊₁.
auto s(index_t i)
BroydenStorage sto
auto s(index_t i) const
BroydenGood(Params params)
length_t current_history() const
Get the number of previous s and y vectors currently stored in the buffer.
bool apply(rvec q)
Apply the inverse Jacobian approximation to the given vector q, i.e.
index_t succ(index_t i) const
Get the next index in the circular buffer of previous s and y vectors.
index_t pred(index_t i) const
Get the previous index in the circular buffer of s and y vectors.
length_t n() const
Get the size of the s and y vectors in the buffer.
const Params & get_params() const
Get the parameters.
length_t history() const
Get the number of previous vectors s and y stored in the buffer.
auto s̃(index_t i)
void foreach_fwd(const F &fun) const
Iterate over the indices in the history buffer, oldest first.
void resize(length_t n)
Re-allocate storage for a problem with a different size.
void reset()
Throw away the approximation and all previous vectors s and y.
bool update_sy(const anymat< VecS > &s, const anymat< VecY > &y, bool forced=false)
Update the inverse Jacobian approximation using the new vectors sₖ = xₖ₊₁ - xₖ and yₖ = pₖ₊₁ - pₖ.
BroydenGood(Params params, length_t n)
Eigen::Ref< const vec > crvec
Default type for immutable references to vectors.
Definition: vec.hpp:18
real_t min_div_abs
Reject update if .
bool restarted
If set to true, the buffer is cleared after memory iterations.
length_t memory
Length of the history to keep.
index_t length_t
Default type for vector sizes.
Definition: vec.hpp:29
bool force_pos_def
If set to true, the inverse Jacobian estimate should remain definite.
Eigen::Index index_t
Default type for vector indices.
Definition: vec.hpp:27
double real_t
Default floating point type.
Definition: vec.hpp:8
Eigen::MatrixBase< Derived > anymat
Generic type for vector and matrix arguments.
Definition: vec.hpp:40
Eigen::Ref< vec > rvec
Default type for mutable references to vectors.
Definition: vec.hpp:16
Parameters for the BroydenGood class.
auto s̃(index_t i) const
auto s(index_t i)
auto s(index_t i) const
void resize(length_t n, length_t history)
Re-allocate storage for a problem with a different size.
length_t n() const
Get the size of the s and s̃ vectors in the buffer.
length_t history() const
Get the number of previous vectors s and s̃ stored in the buffer.
auto s̃(index_t i)
Eigen::Matrix< real_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor > storage_t