|
| | Greenstadt () |
| std::string | name_impl () const |
| Method | method () const |
| void | method (Method t_method) |
| void | enable_one_method () |
| void | enable_two_method () |
| void | set_method (Method t_method) |
| void | update_impl (Vector const &, Vector const &, Matrix const &jacobian_old, Vector const &delta_x_new, Vector const &delta_function_new, Vector const &function_new, Matrix &jacobian_new) |
| | QuasiNewton () |
| std::string | name_impl () const |
| bool | solve_impl (FunctionLambda &&function, JacobianLambda &&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 (FunctionLambda &&function, Vector const &x_ini, Vector &x_sol) |
| | SolverBase () |
| InputType const & | lower_bound () const |
| InputType const & | upper_bound () const |
| void | bounds (InputType const &t_lower_bound, InputType const &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 (FunctionLambda &&function, InputType const &x_ini, InputType &x_sol) |
| bool | rootfind (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, InputType const &x_ini, InputType &x_sol) |
| bool | optimize (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, InputType const &x_ini, InputType &x_sol) |
| std::string | name () const |
template<typename Real,
Integer N>
class Optimist::RootFinder::Greenstadt< Real, N >
Greenstadt's method
Similarly to the Broyden's methods, Greenstadt's methods are quasi-Newton methods that approximates the Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) an update rule. The generic Greenstadt'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.
Based on the update rule, Greenstadt's method is classified into two methods: trivially method 1 and method 2. The update of the Jacobian approximation is performed as
\[ \mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k - \Delta\mathbf{x}_k}{\mathbf{g} \Delta\mathbf{f}_k} \mathbf{g} \text{,}
\]
where \(\Delta\mathbf{x}_k = \mathbf{x}_k - \mathbf{x}_{k-1}\), and \(\Delta\mathbf{f}_k = \mathbf{f}(\mathbf{x}_k) - \mathbf{f}(\mathbf{x}_{k-1})\). The quantity \(\mathbf{g}\) is \(\mathbf{g} = \mathbf{f}_k^\top\) for method 1 and \(\mathbf{g} = \mathbf{H}_{k\top}^{-1} \mathbf{H}_k^{-1} \Delta\mathbf{f}_k^{\top}\) for method 2. The choice of the method is based on the convergence history, switching between the methods to adapt to the problem's behavior.
> For more details on the Greenstadt's methods refer to the reference: Spedicato E., Greenstadt J. On some classes of variationally derived Quasi-Newton methods for systems of nonlinear algebraic equations, Numerische Mathematik, 29 (1978), pp. 363-380.
- Template Parameters
-
| Real | Scalar number type. |
| N | Dimension of the root-finding problem. |