|
| Broyden () |
std::string | name_impl () const |
Method | method () const |
void | method (Method t_method) |
void | enable_good_method () |
void | enable_bad_method () |
void | enable_combined_method () |
void | set_method (Method t_method) |
void | update_impl (Vector const &delta_x_old, Vector const &delta_function_old, Matrix const &jacobian_old, Vector const &delta_x_new, Vector const &delta_function_new, Vector const &, Matrix &jacobian_new) |
| QuasiNewton () |
std::string | name_impl () const |
bool | solve_impl (FunctionWrapper function, JacobianWrapper jacobian, Vector const &x_ini, Vector &x_sol) |
void | update (Vector const &delta_x_old, Vector const &delta_function_old, Matrix const &jacobian_old, Vector const &delta_x_new, Vector const &delta_function_new, Vector const &function_new, Matrix &jacobian_new) |
| RootFinder () |
std::string | name () const |
Integer | jacobian_evaluations () const |
Integer | max_jacobian_evaluations () const |
Integer | hessian_evaluations () const |
Integer | max_hessian_evaluations () const |
bool | solve (FunctionWrapper function, Vector const &x_ini, Vector &x_sol) |
| SolverBase () |
const InputType & | lower_bound () const |
const InputType & | upper_bound () const |
void | bounds (const InputType &t_lower_bound, const InputType &t_upper_bound) |
constexpr Integer | input_dimension () const |
constexpr Integer | output_dimension () const |
Integer | function_evaluations () const |
void | max_function_evaluations (Integer t_max_function_evaluations) |
Integer | iterations () const |
Integer | max_iterations () const |
Real | alpha () const |
Integer | relaxations () const |
Integer | max_relaxations () const |
Real | tolerance () const |
void | verbose_mode (bool t_verbose) |
void | enable_verbose_mode () |
void | disable_verbose_mode () |
void | damped_mode (bool t_damped) |
void | enable_damped_mode () |
void | disable_damped_mode () |
std::string | task () const |
bool | converged () const |
const TraceType & | trace () const |
std::ostream & | ostream () const |
bool | solve (FunctionWrapper function, const InputType &x_ini, InputType &x_sol) |
bool | rootfind (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, const InputType &x_ini, InputType &x_sol) |
bool | optimize (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, const InputType &x_ini, InputType &x_sol) |
std::string | name () const |
template<typename Real,
Integer N>
class Optimist::RootFinder::Broyden< Real, N >
Broyden's method
The Broyden family includes the good, bad, and combined methods. These methods are quasi-Newton approaches that approximate the Jacobian matrix or its inverse. The combined method integrates the updates of the good and bad methods, offering an alternative approach to the classical Newton's method. The Broyden's methods are based on the linearization of the function \(\mathbf{f}(\mathbf{x})\) around the current iterate \(\mathbf{x}_k\), which leads to the linear system
\[ \mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k) \mathbf{h}_k = -\mathbf{f}(\mathbf{x}_k) \text{.}
\]
Here, the Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) is supposed not to be explicitly available or computationally expensive to compute. For this reason, Broyden's methods approximate the Jacobian matrix iteratively. The Broyden's method is defined as
\[ \mathbf{H}_k(\mathbf{x}_k) \mathbf{h}_k = -\mathbf{f}(\mathbf{x}_k) \text{,}
\]
where \(\mathbf{H}_k\) is the (inverse) Jacobian approximation at the \(k\)-th iteration. The generic advancing step is then defined as
\[ \mathbf{x}_{k+1} = \mathbf{x}_k + \alpha_k \mathbf{h}_k \text{.}
\]
where \(\alpha_k\) is a damping coefficient that ensures affine-invariant criteria is satisfied.
What distinguishes Broyden's combined method from the generic Broyden's method is the update of the Jacobian approximation. The Broyden's combined method uses the following update rule
\[ \begin{cases}
\mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k - \Delta\mathbf{x}_k}{\mathbf{C}_g \Delta\mathbf{f}_k} \mathbf{C}_g &
\left\| \displaystyle\frac{\Delta\mathbf{x}_k^\top \Delta\mathbf{x}_{k-1}}{\Delta\mathbf{x}_k^\top \mathbf{H}_k^{-1} \Delta\mathbf{f}_k} \right\| < \left\| \displaystyle\frac{\Delta\mathbf{f}_k^\top \Delta\mathbf{f}_{k-1}}{\Delta\mathbf{f}_k^\top\Delta\mathbf{f}_k} \right\| \\
\mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k-\Delta\mathbf{x}_k}{\mathbf{C}_b \Delta\mathbf{f}_k} \mathbf{C}_b &
\text{otherwise}
\end{cases} \text{,}
\]
with \(\mathbf{C}_g = \Delta\mathbf{x}_k^\top \mathbf{H}_k^{-1}\), \(\mathbf{C}_b = \Delta\mathbf{f}_k^\top\), \(\Delta\mathbf{x}_k = \mathbf{x}_{k+1} - \mathbf{x}_k\), and \(\Delta\mathbf{f}_k = \mathbf{f}(\mathbf{x}_{k+1}) - \mathbf{f}(\mathbf{x}_k)\). Such a rule allows Broyden's combined method to adapt to the problem's behavior, switching between the good and bad methods based on the convergence history.
For more details on the Broyden's methods refer to the references: Broyden C.G., A class of methods for solving nonlinear simultaneous equations, Mathematics of Computation, 19 (1965), pp. 577-593, 10.1090/s0025-5718-1965-0198670-6, and Martinez J., Ochi L., Sobre dois métodos de broyden, Matemática Aplicada e Computacional, IV Congresso Nacional de Matemática Aplicada e Computacional, Rio de Janeiro, Brasil, setembro de 1981.
NOTE: In Optimist, the implemented Broyden's class for the combined method can be used as a good or bad solver by setting the appropriate parameters.
- Template Parameters
-
N | Dimension of the root-finding problem. |
Real | Scalar number type. |