#include <quala/broyden-good.hpp>
Broyden's “Good” method for solving systems of nonlinear equations \( F(x) = 0 \).
\[ \begin{aligned} B_{k+1} &= B_k + \frac{\left(y_k - B_k s_k\right) s_k^\top} {s_k^\top s_k} \\ H_{k+1} &= H_k + \frac{\left(s_k - H_k y_k\right) \left(s_k^\top H_k\right)} {s_k^\top H_k y_k} \\ s_k &\triangleq x_{k+1} - x_k \\ y_k &\triangleq F(x_{k+1}) - F(x_k) \\ \end{aligned} \]
Where \( B_k \) approximates the Jacobian of \( F(x_k) \) and \( H_k \triangleq B_k^{-1} \).
Definition at line 83 of file broyden-good.hpp.
Public Types | |
using | Params = BroydenGoodParams |
Public Member Functions | |
BroydenGood (Params params) | |
BroydenGood (Params params, length_t n) | |
template<class VecS , class VecY > | |
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ₖ. More... | |
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ₖ₊₁. More... | |
bool | apply (rvec q, real_t γ) |
Apply the inverse Jacobian approximation to the given vector q, i.e. More... | |
void | reset () |
Throw away the approximation and all previous vectors s and y. More... | |
void | resize (length_t n) |
Re-allocate storage for a problem with a different size. More... | |
const Params & | get_params () const |
Get the parameters. More... | |
length_t | n () const |
Get the size of the s and y vectors in the buffer. More... | |
length_t | history () const |
Get the number of previous vectors s and y stored in the buffer. More... | |
index_t | succ (index_t i) const |
Get the next index in the circular buffer of previous s and y vectors. More... | |
index_t | pred (index_t i) const |
Get the previous index in the circular buffer of s and y vectors. More... | |
length_t | current_history () const |
Get the number of previous s and y vectors currently stored in the buffer. More... | |
auto | s (index_t i) |
auto | s (index_t i) const |
auto | s̃ (index_t i) |
auto | s̃ (index_t i) const |
template<class F > | |
void | foreach_fwd (const F &fun) const |
Iterate over the indices in the history buffer, oldest first. More... | |
template<class F > | |
void | foreach_rev (const F &fun) const |
Iterate over the indices in the history buffer, newest first. More... | |
Private Attributes | |
BroydenStorage | sto |
index_t | idx = 0 |
bool | full = false |
Params | params |
real_t | latest_γ = NaN |
using Params = BroydenGoodParams |
Definition at line 85 of file broyden-good.hpp.
|
inline |
Definition at line 87 of file broyden-good.hpp.
|
inline |
Update the inverse Jacobian approximation using the new vectors sₖ = xₖ₊₁ - xₖ and yₖ = pₖ₊₁ - pₖ.
Definition at line 177 of file broyden-good.hpp.
Update the inverse Jacobian approximation using the new vectors xₖ₊₁ and pₖ₊₁.
Definition at line 243 of file broyden-good.hpp.
Apply the inverse Jacobian approximation to the given vector q, i.e.
\( q \leftarrow H_k q \). Initial inverse Hessian approximation is set to \( H_0 = \gamma I \). The result is scaled by a factor γ
. If γ
is negative, the result is scaled by \( \frac{s^\top y}{y^\top y} \).
Definition at line 224 of file broyden-good.hpp.
|
inline |
Throw away the approximation and all previous vectors s and y.
Definition at line 164 of file broyden-good.hpp.
|
inline |
Re-allocate storage for a problem with a different size.
Causes a reset.
Definition at line 169 of file broyden-good.hpp.
|
inline |
Get the parameters.
Definition at line 115 of file broyden-good.hpp.
|
inline |
Get the size of the s and y vectors in the buffer.
Definition at line 118 of file broyden-good.hpp.
|
inline |
Get the number of previous vectors s and y stored in the buffer.
Definition at line 120 of file broyden-good.hpp.
Get the next index in the circular buffer of previous s and y vectors.
Definition at line 122 of file broyden-good.hpp.
Get the previous index in the circular buffer of s and y vectors.
Definition at line 124 of file broyden-good.hpp.
|
inline |
Get the number of previous s and y vectors currently stored in the buffer.
Definition at line 127 of file broyden-good.hpp.
|
inline |
Definition at line 129 of file broyden-good.hpp.
|
inline |
|
inline |
Definition at line 131 of file broyden-good.hpp.
|
inline |
|
inline |
Iterate over the indices in the history buffer, oldest first.
Definition at line 136 of file broyden-good.hpp.
|
inline |
Iterate over the indices in the history buffer, newest first.
Definition at line 147 of file broyden-good.hpp.
|
private |
Definition at line 157 of file broyden-good.hpp.
|
private |
Definition at line 158 of file broyden-good.hpp.
|
private |
Definition at line 159 of file broyden-good.hpp.
|
private |
Definition at line 160 of file broyden-good.hpp.
Definition at line 161 of file broyden-good.hpp.