torchmin.least_squares¶
-
torchmin.
least_squares
(fun, x0, bounds=None, method='trf', ftol=1e-08, xtol=1e-08, gtol=1e-08, x_scale=1.0, tr_solver='lsmr', tr_options=None, max_nfev=None, verbose=0)[source]¶ Solve a nonlinear least-squares problem with bounds on the variables.
Given the residual function \(f: \mathcal{R}^n \rightarrow \mathcal{R}^m\), least_squares finds a local minimum of the residual sum-of-squares (RSS) objective:
\[x^* = \underset{x}{\operatorname{arg\,min\,}} \frac{1}{2} ||f(x)||_2^2 \quad \text{subject to} \quad lb \leq x \leq ub\]The solution is found using variants of the Gauss-Newton method, a modification of Newton’s method tailored to RSS problems.
- Parameters
fun (callable) – Function which computes the vector of residuals, with the signature
fun(x)
. The argumentx
passed to this function is a Tensor of shape (n,) (never a scalar, even for n=1). It must allocate and return a 1-D Tensor of shape (m,) or a scalar.x0 (Tensor or float) – Initial guess on independent variables, with shape (n,). If float, it will be treated as a 1-D Tensor with one element.
bounds (2-tuple of Tensor, optional) – Lower and upper bounds on independent variables. Defaults to no bounds. Each Tensor must match the size of x0 or be a scalar, in the latter case a bound will be the same for all variables. Use
inf
with an appropriate sign to disable bounds on all or some variables.method (str, optional) –
Algorithm to perform minimization. Default is ‘trf’.
’trf’ : Trust Region Reflective algorithm, particularly suitable for large sparse problems with bounds. Generally robust method.
’dogbox’ : COMING SOON. dogleg algorithm with rectangular trust regions, typical use case is small problems with bounds. Not recommended for problems with rank-deficient Jacobian.
ftol (float or None, optional) – Tolerance for termination by the change of the cost function. The optimization process is stopped when
dF < ftol * F
, and there was an adequate agreement between a local quadratic model and the true model in the last step. If None, the termination by this condition is disabled. Default is 1e-8.xtol (float or None, optional) – Tolerance for termination by the change of the independent variables. Termination occurs when
norm(dx) < xtol * (xtol + norm(x))
. If None, the termination by this condition is disabled. Default is 1e-8.gtol (float or None, optional) –
Tolerance for termination by the norm of the gradient. Default is 1e-8. The exact condition depends on method used:
For ‘trf’ :
norm(g_scaled, ord=inf) < gtol
, whereg_scaled
is the value of the gradient scaled to account for the presence of the bounds [STIR].For ‘dogbox’ :
norm(g_free, ord=inf) < gtol
, whereg_free
is the gradient with respect to the variables which are not in the optimal state on the boundary.
x_scale (Tensor or 'jac', optional) – Characteristic scale of each variable. Setting x_scale is equivalent to reformulating the problem in scaled variables
xs = x / x_scale
. An alternative view is that the size of a trust region along jth dimension is proportional tox_scale[j]
. Improved convergence may be achieved by setting x_scale such that a step of a given size along any of the scaled variables has a similar effect on the cost function. If set to ‘jac’, the scale is iteratively updated using the inverse norms of the columns of the Jacobian matrix (as described in [JJMore]).max_nfev (None or int, optional) – Maximum number of function evaluations before the termination. Defaults to 100 * n.
tr_solver (str, optional) –
Method for solving trust-region subproblems.
’exact’ is suitable for not very large problems with dense Jacobian matrices. The computational complexity per iteration is comparable to a singular value decomposition of the Jacobian matrix.
’lsmr’ is suitable for problems with sparse and large Jacobian matrices. It uses an iterative procedure for finding a solution of a linear least-squares problem and only requires matrix-vector product evaluations.
tr_options (dict, optional) –
Keyword options passed to trust-region solver.
tr_solver='exact'
: tr_options are ignored.tr_solver='lsmr'
: options for scipy.sparse.linalg.lsmr. Additionally,method='trf'
supports ‘regularize’ option (bool, default is True), which adds a regularization term to the normal equation, which improves convergence if the Jacobian is rank-deficient [Byrd] (eq. 3.4).
verbose (int, optional) –
Level of algorithm’s verbosity.
0 : work silently (default).
1 : display a termination report.
2 : display progress during iterations.
- Returns
result – Result of the optimization routine.
- Return type
OptimizeResult
References
- STIR
M. A. Branch, T. F. Coleman, and Y. Li, “A Subspace, Interior, and Conjugate Gradient Method for Large-Scale Bound-Constrained Minimization Problems,” SIAM Journal on Scientific Computing, Vol. 21, Number 1, pp 1-23, 1999.
- Byrd
R. H. Byrd, R. B. Schnabel and G. A. Shultz, “Approximate solution of the trust region problem by minimization over two-dimensional subspaces”, Math. Programming, 40, pp. 247-263, 1988.
- JJMore
J. J. More, “The Levenberg-Marquardt Algorithm: Implementation and Theory,” Numerical Analysis, ed. G. A. Watson, Lecture Notes in Mathematics 630, Springer Verlag, pp. 105-116, 1977.