|
| RootFinder () |
|
std::string | name () const |
|
Integer | jacobian_evaluations () const |
|
Integer | max_jacobian_evaluations () const |
|
void | max_jacobian_evaluations (Integer t_jacobian_evaluations) |
|
Integer | hessian_evaluations () const |
|
Integer | max_hessian_evaluations () const |
|
void | max_hessian_evaluations (Integer t_hessian_evaluations) |
|
bool | solve (FunctionWrapper function, Vector const &x_ini, Vector &x_sol) |
|
bool | solve (FunctionWrapper function, JacobianWrapper jacobian, Vector const &x_ini, Vector &x_sol) |
|
bool | solve (FunctionWrapper function, JacobianWrapper jacobian, HessianWrapper hessian, Vector const &x_ini, Vector &x_sol) |
|
| Solver () |
|
| Solver (FunctionWrapper function, const InputType &x_ini, InputType &x_sol) |
|
| Solver (FunctionWrapper function, FirstDerivativeWrapper first_derivative, const InputType &x_ini, InputType &x_sol) |
|
| Solver (FunctionWrapper function, FirstDerivativeWrapper first_derivative, SecondDerivativeWrapper second_derivative, const InputType &x_ini, InputType &x_sol) |
|
const InputType & | lower_bound () const |
|
void | lower_bound (const InputType &t_lower_bound) |
|
const InputType & | upper_bound () const |
|
void | upper_bound (const InputType &t_upper_bound) |
|
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 | max_function_evaluations () const |
|
Integer | iterations () const |
|
Integer | max_iterations () const |
|
void | max_iterations (Integer t_max_iterations) |
|
Real | alpha () const |
|
void | alpha (Real t_alpha) |
|
Integer | relaxations () const |
|
Integer | max_relaxations () const |
|
void | max_relaxations (Integer t_max_relaxations) |
|
Real | tolerance () const |
|
void | tolerance (Real t_tolerance) |
|
void | verbose_mode (bool t_verbose) |
|
bool | verbose_mode () const |
|
void | enable_verbose_mode () |
|
void | disable_verbose_mode () |
|
void | damped_mode (bool t_damped) |
|
bool | damped_mode () const |
|
void | enable_damped_mode () |
|
void | disable_damped_mode () |
|
std::string | task () const |
|
void | task (std::string t_task) |
|
bool | converged () const |
|
const TraceType & | trace () const |
|
std::ostream & | ostream () const |
|
void | ostream (std::ostream &t_ostream) |
|
bool | solve (FunctionWrapper function, const InputType &x_ini, InputType &x_sol) |
|
bool | solve (FunctionWrapper function, FirstDerivativeWrapper first_derivative, const InputType &x_ini, InputType &x_sol) |
|
bool | solve (FunctionWrapper function, FirstDerivativeWrapper first_derivative, SecondDerivativeWrapper second_derivative, const InputType &x_ini, InputType &x_sol) |
|
bool | rootfind (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol) |
|
bool | optimize (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol) |
|
std::string | name () const |
|
template<typename Real,
Integer N, typename DerivedSolver>
class Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
Root finder
This section describes the multi-dimensional root-finders implemented in Optimist
. The available root-finding solvers are Newton (derivative) and quasi-Newton (non-derivative) methods. Derivative methods employ the function's derivative to find the root with high accuracy, while non-derivative methods approximate the derivative for improved efficiency in certain scenarios.
Here, the solvers are implemented for solving problems of the form
\[ \mathbf{f}(\mathbf{x}) = \mathbf{0} \text{,} \quad \text{with} \quad \mathbf{f}: \mathbb{R}^{n} \rightarrow \mathbb{R}^{n} \text{,}
\]
which consist in finding the root of the function \(\mathbf{f}\) by iteratively updating the current iterate \(\mathbf{x}_k\) until convergence is achieved. The solvers require the function \(\mathbf{f}\) and its Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) to be provided by the user. The Jacobian matrix can be computed analytically or numerically, depending on the problem's complexity and the user's preference. Alternatively, the Jacobian can be approximated numerically using finite differences, depending on the problem's complexity and the user's preference.
Affine-invariant step
By default, the multi-dimensional root-finders in Optimist
use optionally an affine-invariant advancing step to ensure convergence. The generic advancing step \(\mathbf{h}_k\) is computed as
\[ \mathbf{x}_{k+1} = \mathbf{x}_k + \alpha_k \mathbf{h}_k \text{,}
\]
where \(\alpha_k\) is a damping coefficient that ensures convergence by satisfying
\[ \left\|\mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k)^{-1} \mathbf{f}(\mathbf{x}_{k+1})\right\|
\leq \left(1 - \displaystyle\frac{\alpha_k}{2}\right) \left\|\mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k)^{-1}
\mathbf{f}(\mathbf{x}_k)\right\| = \left(1 - \displaystyle\frac{\alpha_k}{2} \right)
\left\|\mathbf{h}_k\right\| \text{.}
\]
For detailed information on the affine-invariant Newton's method, refer to this link, or the works by P. Deuflhard.
- Template Parameters
-
Real | Scalar number type. |
N | The dimension of the root-finding problem. |
DerivedSolver | Derived solver class. |