naginterfaces.library.opt.lsq_​uncon_​mod_​deriv2_​easy

naginterfaces.library.opt.lsq_uncon_mod_deriv2_easy(m, lsfun2, lshes2, x, data=None, spiked_sorder='C')[source]

lsq_uncon_mod_deriv2_easy is an easy-to-use modified Gauss–Newton algorithm for finding an unconstrained minimum of a sum of squares of nonlinear functions in variables . First and second derivatives are required.

It is intended for functions which are continuous and which have continuous first and second derivatives (although it will usually work even if the derivatives have occasional discontinuities).

For full information please refer to the NAG Library document for e04hy

https://www.nag.com/numeric/nl/nagdoc_29.3/flhtml/e04/e04hyf.html

Parameters
mint

The number of residuals, , and the number of variables, .

lsfun2callable (fvec, fjac) = lsfun2(m, xc, data=None)

You must supply this function to calculate the vector of values and the Jacobian matrix of first derivatives at any point .

It should be tested separately before being used in conjunction with lsq_uncon_mod_deriv2_easy (see the E04 Introduction).

Parameters
mint

, the numbers of residuals.

xcfloat, ndarray, shape

The point at which the values of the and the are required.

dataarbitrary, optional, modifiable in place

User-communication data for callback functions.

Returns
fvecfloat, array-like, shape

must be set to the value of at the point , for .

fjacfloat, array-like, shape

must be set to the value of at the point , for , for .

lshes2callable b = lshes2(fvec, xc, lb, data=None)

You must supply this function to calculate the elements of the symmetric matrix

at any point , where is the Hessian matrix of .

It should be tested separately before being used in conjunction with lsq_uncon_mod_deriv2_easy (see the E04 Introduction).

Parameters
fvecfloat, ndarray, shape

The value of the residual at the point , for , so that the values of the can be used in the calculation of the elements of .

xcfloat, ndarray, shape

The point at which the elements of are to be evaluated.

lbint

The length of the array .

dataarbitrary, optional, modifiable in place

User-communication data for callback functions.

Returns
bfloat, array-like, shape

Must contain the lower triangle of the matrix , evaluated at the point , stored by rows. (The upper triangle is not required because the matrix is symmetric.) More precisely, must contain evaluated at the point , for , for .

xfloat, array-like, shape

must be set to a guess at the th component of the position of the minimum, for . The function checks and at the starting point and so is more likely to detect any error in your functions if the initial are nonzero and mutually distinct.

dataarbitrary, optional

User-communication data for callback functions.

spiked_sorderstr, optional

If in is spiked (i.e., has unit extent in all but one dimension, or has size ), selects the storage order to associate with it in the NAG Engine:

spiked_sorder =

row-major storage will be used;

spiked_sorder =

column-major storage will be used.

Returns
xfloat, ndarray, shape

The lowest point found during the calculations. Thus, if no exception or warning is raised on exit, is the th component of the position of the minimum.

fsumsqfloat

The value of the sum of squares, , corresponding to the final point stored in .

commdict, communication object

Communication structure.

Raises
NagValueError
(errno )

On entry, and .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

There have been calls to .

(errno )

Failure in computing SVD of Jacobian matrix.

(errno )

It is very likely that you have made an error in forming the derivatives in .

(errno )

It is very likely that you have made an error in setting up the array in .

Warns
NagAlgorithmicWarning
(errno )

The conditions for a minimum have not all been satisfied, but a lower point could not be found.

(errno )

It is probable that a local minimum has been found, but it cannot be guaranteed.

(errno )

It is possible that a local minimum has been found, but it cannot be guaranteed.

(errno )

It is unlikely that a local minimum has been found.

(errno )

It is very unlikely that a local minimum has been found.

Notes

No equivalent traditional C interface for this routine exists in the NAG Library.

lsq_uncon_mod_deriv2_easy is similar to the function LSSDN2 in the NPL Algorithms Library. It is applicable to problems of the form:

where and . (The functions are often referred to as ‘residuals’.)

You must supply a function to evaluate the residuals and their first derivatives at any point , and a function to evaluate the elements of the second derivative term of the Hessian matrix of .

Before attempting to minimize the sum of squares, the algorithm checks the functions for consistency. Then, from a starting point supplied by you, a sequence of points is generated which is intended to converge to a local minimum of the sum of squares. These points are generated using estimates of the curvature of .

References

Gill, P E and Murray, W, 1978, Algorithms for the solution of the nonlinear least squares problem, SIAM J. Numer. Anal. (15), 977–992