quala main
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 /// Powell's trick, damping, prevents nonsingularity.
58};
59
60/**
61 * Broyden's “Good” method for solving systems of nonlinear equations
62 * @f$ F(x) = 0 @f$.
63 *
64 * @f[ \begin{aligned}
65 * B_{k+1} &= B_k + \frac{\left(y_k - B_k s_k\right) s_k^\top}
66 * {s_k^\top s_k} \\
67 * H_{k+1} &= H_k + \frac{\left(s_k - H_k y_k\right)
68 * \left(s_k^\top H_k\right)}
69 * {s_k^\top H_k y_k} \\
70 * s_k &\triangleq x_{k+1} - x_k \\
71 * y_k &\triangleq F(x_{k+1}) - F(x_k) \\
72 * \end{aligned} @f]
73 * Where @f$ B_k @f$ approximates the Jacobian of @f$ F(x_k) @f$ and
74 * @f$ H_k \triangleq B_k^{-1} @f$.
75 *
76 * @todo Damping.
77 *
78 * @ingroup accelerators-grp
79 */
81 public:
83
86
87 /// Update the inverse Jacobian approximation using the new vectors
88 /// sₖ = xₖ₊₁ - xₖ and yₖ = pₖ₊₁ - pₖ.
89 template <class VecS, class VecY>
90 bool update_sy(const anymat<VecS> &s, const anymat<VecY> &y,
91 bool forced = false);
92
93 /// Update the inverse Jacobian approximation using the new vectors xₖ₊₁
94 /// and pₖ₊₁.
95 bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁,
96 bool forced = false);
97
98 /// Apply the inverse Jacobian approximation to the given vector q, i.e.
99 /// @f$ q \leftarrow H_k q @f$.
100 /// Initial inverse Hessian approximation is set to @f$ H_0 = \gamma I @f$.
101 /// The result is scaled by a factor @p γ. If @p γ is negative, the result
102 /// is scaled by @f$ \frac{s^\top y}{y^\top y} @f$.
103 bool apply(rvec q, real_t γ);
104
105 /// Throw away the approximation and all previous vectors s and y.
106 void reset();
107 /// Re-allocate storage for a problem with a different size. Causes
108 /// a @ref reset.
109 void resize(length_t n);
110
111 /// Get the parameters.
112 const Params &get_params() const { return params; }
113
114 /// Get the size of the s and y vectors in the buffer.
115 length_t n() const { return sto.n(); }
116 /// Get the number of previous vectors s and y stored in the buffer.
117 length_t history() const { return sto.history(); }
118 /// Get the next index in the circular buffer of previous s and y vectors.
119 index_t succ(index_t i) const { return i + 1 < history() ? i + 1 : 0; }
120 /// Get the previous index in the circular buffer of s and y vectors.
121 index_t pred(index_t i) const { return i > 0 ? i - 1 : history() - 1; }
122 /// Get the number of previous s and y vectors currently stored in the
123 /// buffer.
124 length_t current_history() const { return full ? history() : idx; }
125
126 auto s(index_t i) { return sto.s(i); }
127 auto s(index_t i) const { return sto.s(i); }
128 auto (index_t i) { return sto.(i); }
129 auto (index_t i) const { return sto.(i); }
130
131 /// Iterate over the indices in the history buffer, oldest first.
132 template <class F>
133 void foreach_fwd(const F &fun) const {
134 if (full)
135 for (index_t i = idx; i < history(); ++i)
136 fun(i);
137 if (idx)
138 for (index_t i = 0; i < idx; ++i)
139 fun(i);
140 }
141
142 /// Iterate over the indices in the history buffer, newest first.
143 template <class F>
144 void foreach_rev(const F &fun) const {
145 if (idx)
146 for (index_t i = idx; i-- > 0;)
147 fun(i);
148 if (full)
149 for (index_t i = history(); i-- > idx;)
150 fun(i);
151 }
152
153 private:
156 bool full = false;
159};
160
161inline void BroydenGood::reset() {
162 idx = 0;
163 full = false;
164}
165
167 if (params.memory < 1)
168 throw std::invalid_argument("BroydenGood::Params::memory must be >= 1");
170 reset();
171}
172
173template <class VecS, class VecY>
175 bool forced) {
176 // Restart if the buffer is full
177 if (full && params.restarted) {
178 full = false;
179 assert(idx == 0);
180 }
181
182 auto &&r = sto.work();
183 // Compute r = r₍ₘ₋₁₎ = Hₖ yₖ
184 r = yₖ; // r₍₋₁₎ = yₖ
185 foreach_fwd([&](index_t i) {
186 r += (i) * r.dot(s(i)); // r₍ᵢ₎ = r₍ᵢ₋₁₎ + s̃₍ᵢ₎〈r₍ᵢ₋₁₎, s₍ᵢ₎〉
187 });
188 const real_t sᵀHy = sₖ.dot(r);
189 const real_t a_sᵀHy = params.force_pos_def ? sᵀHy : std::abs(sᵀHy);
190
191 // Check if update is accepted
192 if (!forced && a_sᵀHy < params.min_div_abs)
193 return false;
194
195 // Compute damping
196 real_t damp = 1 / sᵀHy;
198 real_t γ = sᵀHy / sₖ.squaredNorm();
199 real_t sgn_γ = γ >= 0 ? 1 : -1;
200 real_t γθ̅ = params.force_pos_def ? γ * θ̅ : sgn_γ * θ̅;
201 real_t a_γ = params.force_pos_def ? γ : std::abs(γ);
202 real_t θ = a_γ >= θ̅ ? 1 // no damping
203 : (1 - γθ̅) / (1 - γ);
204 damp = θ / (sₖ.squaredNorm() * (1 - θ + θ * γ));
205 }
206
207 // Store the new vectors
208 sto.s(idx) = sₖ;
209 sto.(idx) = damp * (sₖ - r);
210 latest_γ = sₖ.dot(yₖ) / yₖ.squaredNorm();
211
212 // Increment the index in the circular buffer
213 idx = succ(idx);
214 full |= idx == 0;
215
216 return true;
217}
218
219inline bool BroydenGood::apply(rvec q, real_t γ) {
220 // Only apply if we have previous vectors s and y
221 if (idx == 0 && not full)
222 return false;
223
224 if (γ < 0)
225 γ = latest_γ;
226 if (γ != 1)
227 q *= γ;
228
229 // Compute q = q₍ₘ₋₁₎ = Hₖ q
230 // q₍₋₁₎ = q
231 foreach_fwd([&](index_t i) {
232 q += (i) * q.dot(s(i)); // q₍ᵢ₎ = q₍ᵢ₋₁₎ + s̃₍ᵢ₎〈q₍ᵢ₋₁₎, s₍ᵢ₎〉
233 });
234
235 return true;
236}
237
238inline bool BroydenGood::update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁,
239 bool forced) {
240 const auto s = xₖ₊₁ - xₖ;
241 const auto y = pₖ₊₁ - pₖ;
242 return update_sy(s, y, forced);
243}
244
245} // 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.
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)
bool apply(rvec q, real_t γ)
Apply the inverse Jacobian approximation to the given vector q, i.e.
real_t powell_damping_factor
Powell's trick, damping, prevents nonsingularity.
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.
constexpr real_t NaN
Not a number.
Definition: vec.hpp:45
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