NAG FL Interface
e04gbf (lsq_​uncon_​quasi_​deriv_​comp)

Note: this routine is deprecated. Replaced by e04ggf.
e04ggf is part of the new NAG optimization modelling suite (see Section 3.1 in the E04 Chapter Introduction), therefore the definition of the nonlinear residual function values and gradients need to be split into two separate subroutines. e04ggf offers a significant improvement in performance over e04gbf as well as additional functionality, such as the addition of variable bounds and user-evaluation recovery, amongst many others.
Callbacks

Old Code

       Subroutine lsqfun(iflag,m,n,xc,fvec,fjac,ldfjac,iw,liw,w,lw)

!        .. Implicit None Statement ..
         Implicit None
!        .. Scalar Arguments ..
         Integer, Intent (Inout)        :: iflag
         Integer, Intent (In)           :: ldfjac, liw, lw, m, n
!        .. Array Arguments ..
         Real (Kind=nag_wp), Intent (Inout) :: fjac(ldfjac,n), w(lw)
         Real (Kind=nag_wp), Intent (Out) :: fvec(m)
         Real (Kind=nag_wp), Intent (In) :: xc(n)
         Integer, Intent (Inout)        :: iw(liw)
!        .. Executable Statements ..

!        Compute residuals
         fvec(:) = ...

         If (iflag/=0) Then
!          Compute Jacobian of residuals
           fjac(:,:) = ...
         End If

         Return
       End Subroutine lsqfun

New Code

       Subroutine lsqfun(nvar,x,nres,rx,inform,iuser,ruser,cpuser)

!        .. Use Statements ..
         Use, Intrinsic                 :: iso_c_binding, Only: c_ptr
!        .. Scalar Arguments ..
         Type (c_ptr), Intent (In)      :: cpuser
         Integer, Intent (Inout)        :: inform
         Integer, Intent (In)           :: nres, nvar
!        .. Array Arguments ..
         Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
         Real (Kind=nag_wp), Intent (Out) :: rx(nres)
         Real (Kind=nag_wp), Intent (In) :: x(nvar)
         Integer, Intent (Inout)        :: iuser(*)
!        .. Executable Statements ..

!        Compute residuals
         rx(:) = ...

!        Inform evaluation was successful
         inform = 0
         Return
       End Subroutine lsqfun

       Subroutine lsqgrd(nvar,x,nres,nnzrd,rdx,inform,iuser,ruser,cpuser)

!        .. Use Statements ..
         Use, Intrinsic                 :: iso_c_binding, Only: c_ptr
!        .. Scalar Arguments ..
         Type (c_ptr), Intent (In)      :: cpuser
         Integer, Intent (Inout)        :: inform
         Integer, Intent (In)           :: nnzrd, nres, nvar
!        .. Array Arguments ..
         Real (Kind=nag_wp), Intent (Inout) :: rdx(nnzrd), ruser(*)
         Real (Kind=nag_wp), Intent (In) :: x(nvar)
         Integer, Intent (Inout)        :: iuser(*)
!        .. Executable Statements ..

!        Compute gradient of residuals
         rdx(:) = ...

!        Inform evaluation was successful
         inform = 0
         Return
       End Subroutine lsqgrd
Main Call

Old Code

       ifail = -1
       Call e04gbf(m,n,e04hev,lsqfun,e04fdz,iprint,maxcal,eta,xtol,stepmx,x,    &
         fsumsq,fvec,fjac,ldfjac,s,v,ldv,niter,nf,iw,liw,w,lw,ifail)

New Code

!      .. Use Statements ..
       Use, Intrinsic                   :: iso_c_binding, Only: c_null_ptr, c_ptr

!      ...

!      Initialize problem with n variables
       ifail = 0
       Call e04raf(handle,n,ifail)

!      Define nonlinear least-square problem with m residuals and dense structure
       isparse = 0
       nnzrd = 0
       Call e04rmf(handle,m,isparse,nnzrd,irowrd,icolrd,ifail)

!      Solve the problem
       cpuser = c_null_ptr
       ifail = -1
       Call e04ggf(handle,lsqfun,lsqgrd,e04ggu,e04ggv,e04ffu,n,x,m,rx,rinfo,    &
         stats,iuser,ruser,cpuser,ifail)

       iter = Int(stats(1))
       lsqf = rinfo(1)
       error = sqrt(2.0_nag_wp*lsqf)

!      Free the handle memory
       ifail = 0
       Call e04rzf(handle,ifail)
Settings help

FL Name Style:


FL Specification Language:


1 Purpose

e04gbf is a comprehensive quasi-Newton algorithm for finding an unconstrained minimum of a sum of squares of m nonlinear functions in n variables (mn). First derivatives are required. The routine is intended for functions which have continuous first and second derivatives (although it will usually work even if the derivatives have occasional discontinuities).

2 Specification

Fortran Interface
Subroutine e04gbf ( m, n, lsqlin, lsqfun, lsqmon, iprint, maxcal, eta, xtol, stepmx, x, fsumsq, fvec, fjac, ldfjac, s, v, ldv, niter, nf, iw, liw, w, lw, ifail)
Integer, Intent (In) :: m, n, iprint, maxcal, ldfjac, ldv, liw, lw
Integer, Intent (Inout) :: iw(liw), ifail
Integer, Intent (Out) :: niter, nf
Real (Kind=nag_wp), Intent (In) :: eta, xtol, stepmx
Real (Kind=nag_wp), Intent (Inout) :: x(n), fjac(ldfjac,n), v(ldv,n), w(lw)
Real (Kind=nag_wp), Intent (Out) :: fsumsq, fvec(m), s(n)
External :: lsqlin, lsqfun, lsqmon
C Header Interface
#include <nag.h>
void  e04gbf_ (const Integer *m, const Integer *n,
void (NAG_CALL *lsqlin)(Integer *selct),
void (NAG_CALL *lsqfun)(Integer *iflag, const Integer *m, const Integer *n, const double xc[], double fvec[], double fjac[], const Integer *ldfjac, Integer iw[], const Integer *liw, double w[], const Integer *lw),
void (NAG_CALL *lsqmon)(const Integer *m, const Integer *n, const double xc[], const double fvec[], const double fjac[], const Integer *ldfjac, const double s[], const Integer *igrade, const Integer *niter, const Integer *nf, Integer iw[], const Integer *liw, double w[], const Integer *lw),
const Integer *iprint, const Integer *maxcal, const double *eta, const double *xtol, const double *stepmx, double x[], double *fsumsq, double fvec[], double fjac[], const Integer *ldfjac, double s[], double v[], const Integer *ldv, Integer *niter, Integer *nf, Integer iw[], const Integer *liw, double w[], const Integer *lw, Integer *ifail)
The routine may be called by the names e04gbf or nagf_opt_lsq_uncon_quasi_deriv_comp.

3 Description

e04gbf is essentially identical to the subroutine LSQFDQ in the NPL Algorithms Library. It is applicable to problems of the form:
MinimizeF(x)=i=1m[fi(x)]2  
where x = (x1,x2,,xn) T and mn. (The functions fi(x) are often referred to as ‘residuals’.)
You must supply a subroutine to calculate the values of the fi(x) and their first derivatives fi xj at any point x.
From a starting point x (1) supplied by you, the routine generates a sequence of points x (2) ,x (3) ,, which is intended to converge to a local minimum of F(x). The sequence of points is given by
x (k+1) =x (k) +α(k)p (k)  
where the vector p (k) is a direction of search, and α(k) is chosen such that F(x (k) +α(k)p (k) ) is approximately a minimum with respect to α(k).
The vector p (k) used depends upon the reduction in the sum of squares obtained during the last iteration. If the sum of squares was sufficiently reduced, then p (k) is the Gauss–Newton direction; otherwise the second derivatives of the fi(x) are taken into account using a quasi-Newton updating scheme.
The method is designed to ensure that steady progress is made whatever the starting point, and to have the rapid ultimate convergence of Newton's method.

4 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

5 Arguments

1: m Integer Input
2: n Integer Input
On entry: the number m of residuals, fi(x), and the number n of variables, xj.
Constraint: 1nm.
3: lsqlin Subroutine, supplied by the NAG Library. External Procedure
lsqlin enables you to specify whether the linear minimizations (i.e., minimizations of F(x (k) +α (k) p (k) ) with respect to α(k)) are to be performed by a routine which just requires the evaluation of the fi(x) (e04fcv), or by a routine which also requires the first derivatives of the fi(x) (e04hev).
It will often be possible to evaluate the first derivatives of the residuals in about the same amount of computer time that is required for the evaluation of the residuals themselves – if this is so, then e04gbf should be called with routine e04hev as the argument lsqlin. However, if the evaluation of the derivatives takes more than about four times as long as the evaluation of the residuals, then e04fcv will usually be preferable. If in doubt, use e04hev as it is slightly more robust.
Whichever subroutine is used, must be declared as EXTERNAL in the subroutine from which e04gbf is called.
4: lsqfun Subroutine, supplied by the user. External Procedure
lsqfun must calculate the vector of values fi(x) and Jacobian matrix of first derivatives fi xj at any point x. (However, if you do not wish to calculate the residuals or first derivatives at a particular x, there is the option of setting an argument to cause e04gbf to terminate immediately.)
The specification of lsqfun is:
Fortran Interface
Subroutine lsqfun ( iflag, m, n, xc, fvec, fjac, ldfjac, iw, liw, w, lw)
Integer, Intent (In) :: m, n, ldfjac, liw, lw
Integer, Intent (Inout) :: iflag, iw(liw)
Real (Kind=nag_wp), Intent (In) :: xc(n)
Real (Kind=nag_wp), Intent (Inout) :: fjac(ldfjac,n), w(lw)
Real (Kind=nag_wp), Intent (Out) :: fvec(m)
C Header Interface
void  lsqfun (Integer *iflag, const Integer *m, const Integer *n, const double xc[], double fvec[], double fjac[], const Integer *ldfjac, Integer iw[], const Integer *liw, double w[], const Integer *lw)
Important: the dimension declaration for fjac must contain the variable ldfjac, not an integer constant.
1: iflag Integer Input/Output
On entry: will be set to 0, 1 or 2.
iflag=0
Indicates that only the residuals need to be evaluated
iflag=1
Indicates that only the Jacobian matrix needs to be evaluated
iflag=2
Indicates that both the residuals and the Jacobian matrix must be calculated.
If e04hev is used as e04gbf's lsqlin, lsqfun will always be called with iflag set to 2.
On exit: if it is not possible to evaluate the fi(x) or their first derivatives at the point given in xc (or if it is wished to stop the calculations for any other reason), you should reset iflag to some negative number and return control to e04gbf. e04gbf will then terminate immediately, with ifail set to your setting of iflag.
2: m Integer Input
On entry: m, the number of residuals.
3: n Integer Input
On entry: n, the number of variables.
4: xc(n) Real (Kind=nag_wp) array Input
On entry: the point x at which the values of the fi and the fi xj are required.
5: fvec(m) Real (Kind=nag_wp) array Output
On exit: unless iflag=1 on entry, or iflag is reset to a negative number, fvec(i) must contain the value of fi at the point x, for i=1,2,,m.
6: fjac(ldfjac,n) Real (Kind=nag_wp) array Output
On exit: unless iflag=0 on entry, or iflag is reset to a negative number, fjac(i,j) must contain the value of fi xj at the point x, for i=1,2,,m and j=1,2,,n.
7: ldfjac Integer Input
On entry: the first dimension of the array fjac, set to m by e04gbf.
8: iw(liw) Integer array Workspace
9: liw Integer Input
10: w(lw) Real (Kind=nag_wp) array Workspace
11: lw Integer Input
lsqfun is called with e04gbf's arguments iw, liw, w, lw as these arguments. They are present so that, when other library routines require the solution of a minimization subproblem, constants needed for the evaluation of residuals can be passed through iw and w. Similarly, you could pass quantities to lsqfun from the segment which calls e04gbf by using partitions of iw and w beyond those used as workspace by e04gbf. However, because of the danger of mistakes in partitioning, it is recommended that you should pass information to lsqfun via COMMON global variables and not use iw or w at all. In any case you must not change the elements of iw and w used as workspace by e04gbf.
lsqfun must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which e04gbf is called. Arguments denoted as Input must not be changed by this procedure.
Note: lsqfun should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by e04gbf. If your code inadvertently does return any NaNs or infinities, e04gbf is likely to produce unexpected results.
lsqfun should be tested separately before being used in conjunction with e04gbf.
5: lsqmon Subroutine, supplied by the NAG Library or the user. External Procedure
If iprint0, you must supply lsqmon which is suitable for monitoring the minimization process. lsqmon must not change the values of any of its arguments.
If iprint<0, the NAG Library dummy routine e04fdz can be used as lsqmon.
The specification of lsqmon is:
Fortran Interface
Subroutine lsqmon ( m, n, xc, fvec, fjac, ldfjac, s, igrade, niter, nf, iw, liw, w, lw)
Integer, Intent (In) :: m, n, ldfjac, igrade, niter, nf, liw, lw
Integer, Intent (Inout) :: iw(liw)
Real (Kind=nag_wp), Intent (In) :: xc(n), fvec(m), fjac(ldfjac,n), s(n)
Real (Kind=nag_wp), Intent (Inout) :: w(lw)
C Header Interface
void  lsqmon (const Integer *m, const Integer *n, const double xc[], const double fvec[], const double fjac[], const Integer *ldfjac, const double s[], const Integer *igrade, const Integer *niter, const Integer *nf, Integer iw[], const Integer *liw, double w[], const Integer *lw)
Important: the dimension declaration for fjac must contain the variable ldfjac, not an integer constant.
1: m Integer Input
On entry: m, the numbers of residuals.
2: n Integer Input
On entry: n, the numbers of variables.
3: xc(n) Real (Kind=nag_wp) array Input
On entry: the coordinates of the current point x.
4: fvec(m) Real (Kind=nag_wp) array Input
On entry: the values of the residuals fi at the current point x.
5: fjac(ldfjac,n) Real (Kind=nag_wp) array Input
On entry: fjac(i,j) contains the value of fi xj at the current point x, for i=1,2,,m and j=1,2,,n.
6: ldfjac Integer Input
On entry: the first dimension of the array fjac as declared in the (sub)program from which e04gbf is called.
7: s(n) Real (Kind=nag_wp) array Input
On entry: the singular values of the current Jacobian matrix. Thus s may be useful as information about the structure of your problem.
8: igrade Integer Input
On entry: e04gbf estimates the dimension of the subspace for which the Jacobian matrix can be used as a valid approximation to the curvature (see Gill and Murray (1978)). This estimate is called the grade of the Jacobian matrix, and igrade gives its current value.
9: niter Integer Input
On entry: the number of iterations which have been performed in e04gbf.
10: nf Integer Input
On entry: the number of evaluations of the residuals. (If e04hev is used as lsqlin, nf is also the number of evaluations of the Jacobian matrix.)
11: iw(liw) Integer array Workspace
12: liw Integer Input
13: w(lw) Real (Kind=nag_wp) array Workspace
14: lw Integer Input
As in lsqfun, these arguments correspond to the arguments iw, liw, w, lw of e04gbf. They are included in lsqmon's argument list primarily for when e04gbf is called by other library routines.
lsqmon must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which e04gbf is called. Arguments denoted as Input must not be changed by this procedure.
Note:  you should normally print the sum of squares of residuals, so as to be able to examine the sequence of values of F(x) mentioned in Section 7. It is usually helpful to also print xc, the gradient of the sum of squares, niter and nf.
6: iprint Integer Input
On entry: the frequency with which lsqmon is to be called.
iprint>0
lsqmon is called once every iprint iterations and just before exit from e04gbf.
iprint=0
lsqmon is just called at the final point.
iprint<0
lsqmon is not called at all.
iprint should normally be set to a small positive number.
Suggested value: iprint=1.
7: maxcal Integer Input
On entry: enables you to limit the number of times that lsqfun is called by e04gbf. There will be an error exit (see Section 6) after maxcal calls of lsqfun.
Suggested values:
  • maxcal=75×n if e04fcv is used as lsqlin,
  • maxcal=50×n if e04hev is used as lsqlin.
Constraint: maxcal1.
8: eta Real (Kind=nag_wp) Input
On entry: every iteration of e04gbf involves a linear minimization (i.e., minimization of F(x (k) +α (k) p (k) ) with respect to α (k) ). eta specifies how accurately these linear minimizations are to be performed. The minimum with respect to α (k) will be located more accurately for small values of eta (say, 0.01) than for large values (say, 0.9).
Although accurate linear minimizations will generally reduce the number of iterations performed by e04gbf, they will increase the number of calls of lsqfun made every iteration. On balance it is usually more efficient to perform a low accuracy minimization.
Suggested values:
  • eta=0.9 if n>1 and e04hev is used as lsqlin,
  • eta=0.5 if n>1 and e04fcv is uses as lsqlin,
  • eta=0.0 if n=1.
Constraint: 0.0eta<1.0.
9: xtol Real (Kind=nag_wp) Input
On entry: the accuracy in x to which the solution is required.
If xtrue is the true value of x at the minimum, then xsol, the estimated position before a normal exit, is such that
xsol-xtrue<xtol×(1.0+xtrue),  
where y=j=1nyj2. For example, if the elements of xsol are not much larger than 1.0 in modulus and if xtol=1.0E−5, then xsol is usually accurate to about five decimal places. (For further details see Section 7.)
If F(x) and the variables are scaled roughly as described in Section 9 and ε is the machine precision, then a setting of order xtol=ε will usually be appropriate. If xtol is set to 0.0 or some positive value less than 10ε, e04gbf will use 10ε instead of xtol, since 10ε is probably the smallest reasonable setting.
Constraint: xtol0.0.
10: stepmx Real (Kind=nag_wp) Input
On entry: an estimate of the Euclidean distance between the solution and the starting point supplied by you. (For maximum efficiency, a slight overestimate is preferable.)
e04gbf will ensure that, for each iteration,
j=1n(xj (k) -xj (k-1) )2 (stepmx) 2  
where k is the iteration number. Thus, if the problem has more than one solution, e04gbf is most likely to find the one nearest to the starting point. On difficult problems, a realistic choice can prevent the sequence of x (k) entering a region where the problem is ill-behaved and can help avoid overflow in the evaluation of F(x). However, an underestimate of stepmx can lead to inefficiency.
Suggested value: stepmx=100000.0.
Constraint: stepmxxtol.
11: x(n) Real (Kind=nag_wp) array Input/Output
On entry: x(j) must be set to a guess at the jth component of the position of the minimum, for j=1,2,,n.
On exit: the final point x (k) . Thus, if ifail=0 on exit, x(j) is the jth component of the estimated position of the minimum.
12: fsumsq Real (Kind=nag_wp) Output
On exit: the value of F(x), the sum of squares of the residuals fi(x), at the final point given in x.
13: fvec(m) Real (Kind=nag_wp) array Output
On exit: the value of the residual fi(x) at the final point given in x, for i=1,2,,m.
14: fjac(ldfjac,n) Real (Kind=nag_wp) array Output
On exit: the value of the first derivative fi xj evaluated at the final point given in x, for i=1,2,,m and j=1,2,,n.
15: ldfjac Integer Input
On entry: the first dimension of the array fjac as declared in the (sub)program from which e04gbf is called.
Constraint: ldfjacm.
16: s(n) Real (Kind=nag_wp) array Output
On exit: the singular values of the Jacobian matrix at the final point. Thus s may be useful as information about the structure of your problem.
17: v(ldv,n) Real (Kind=nag_wp) array Output
On exit: the matrix V associated with the singular value decomposition
J=USVT  
of the Jacobian matrix at the final point, stored by columns. This matrix may be useful for statistical purposes, since it is the matrix of orthonormalized eigenvectors of JTJ.
18: ldv Integer Input
On entry: the first dimension of the array v as declared in the (sub)program from which e04gbf is called.
Constraint: ldvn.
19: niter Integer Output
On exit: the number of iterations which have been performed in e04gbf.
20: nf Integer Output
On exit: the number of times that the residuals have been evaluated (i.e., the number of calls of lsqfun). If e04hev is used as lsqlin, nf is also the number of times that the Jacobian matrix has been evaluated.
21: iw(liw) Integer array Workspace
22: liw Integer Input
On entry: the dimension of the array iw as declared in the (sub)program from which e04gbf is called.
Constraint: liw1.
23: w(lw) Real (Kind=nag_wp) array Workspace
24: lw Integer Input
On entry: the dimension of the array w as declared in the (sub)program from which e04gbf is called.
Constraints:
  • if n>1, lw7×n+m×n+2×m+n×n;
  • if n=1, lw9+3×m.
25: ifail Integer Input/Output
On entry: ifail must be set to 0, −1 or 1 to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of 0 causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of −1 means that an error message is printed while a value of 1 means that it is not.
If halting is not appropriate, the value −1 or 1 is recommended. If message printing is undesirable, then the value 1 is recommended. Otherwise, the value −1 is recommended since useful values can be provided in some output arguments even when ifail0 on exit. When the value -1 or 1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6 Error Indicators and Warnings

If on entry ifail=0 or −1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
Note: in some cases e04gbf may return useful information.
ifail=1
On entry, eta=value.
Constraint: 0.0eta<1.0.
On entry, ldfjac=value and m=value.
Constraint: ldfjacm.
On entry, ldv=value and n=value.
Constraint: ldvn.
On entry, liw=value.
Constraint: liw1.
On entry, m=value and n=value.
Constraint: mn.
On entry, maxcal=value.
Constraint: maxcal1.
On entry, n=value.
Constraint: n1.
On entry, n=1 and lw=value.
Constraint: if n=1 then lw9+3×m; that is, value.
On entry, n>1 and lw=value.
Constraint: if n>1 then lw7×n+m×n+2×m+n2; that is, value.
On entry, stepmx=value and xtol=value.
Constraint: stepmxxtol.
On entry, xtol=value.
Constraint: xtol0.0.
ifail=2
There have been maxcal=value calls to lsqfun.
If steady reductions in the sum of squares, F(x), were monitored up to the point where this exit occurred, then the exit probably occurred simply because maxcal was set too small, so the calculations should be restarted from the final point held in x. This exit may also indicate that F(x) has no minimum.
The error may also be caused by mistakes in lsqfun, by the formulation of the problem or by an awkward function. If there are no such mistakes, it is worth restarting the calculations from a different starting point (not the point at which the failure occurred) in order to avoid the region which caused the failure.
ifail=3
The conditions for a minimum have not all been satisfied, but a lower point could not be found. See Section 7 for further information.
The error may also be caused by mistakes in lsqfun, by the formulation of the problem or by an awkward function. If there are no such mistakes, it is worth restarting the calculations from a different starting point (not the point at which the failure occurred) in order to avoid the region which caused the failure.
ifail=4
Failure in computing SVD of Jacobian matrix.
It may be worth applying e04gbf again starting with an initial approximation which is not too close to the point at which the failure occurred.
The error may also be caused by mistakes in lsqfun, by the formulation of the problem or by an awkward function. If there are no such mistakes, it is worth restarting the calculations from a different starting point (not the point at which the failure occurred) in order to avoid the region which caused the failure.
ifail<0
User requested termination by setting iflag negative in lsqfun.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.

7 Accuracy

A successful exit (ifail=0) is made from e04gbf when (B1, B2 and B3) or B4 or B5 hold, where
B1 α (k) ×p (k) <(xtol+ε)×(1.0+x (k) ) B2 |F (k) -F (k-1) |< (xtol+ε) 2×(1.0+F (k) ) B3 g (k) <ε1/3×(1.0+F (k) ) B4 F (k) <ε2 B5 g (k) < (ε×F (k) ) 1/2  
and where . and ε are as defined in xtol, and F (k) and g (k) are the values of F(x) and its vector of first derivatives at x (k) .
If ifail=0, then the vector in x on exit, xsol, is almost certainly an estimate of xtrue, the position of the minimum to the accuracy specified by xtol.
If ifail=3, then xsol may still be a good estimate of xtrue, but to verify this you should make the following checks. If
  1. (a)the sequence {F(x (k) )} converges to F(xsol) at a superlinear or a fast linear rate, and
  2. (b) g(xsol)T g(xsol) < 10 ε where T denotes transpose, then it is almost certain that xsol is a close approximation to the minimum.
When (b) is true, then usually F(xsol) is a close approximation to F(xtrue). The values of F(x (k) ) can be calculated in lsqmon, and the vector g(xsol) can be calculated from the contents of fvec and fjac on exit from e04gbf.
Further suggestions about confirmation of a computed solution are given in the E04 Chapter Introduction.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
e04gbf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
e04gbf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

The number of iterations required depends on the number of variables, the number of residuals, the behaviour of F(x), the accuracy demanded and the distance of the starting point from the solution. The number of multiplications performed per iteration of e04gbf varies, but for mn is approximately n×m2+O(n3). In addition, each iteration makes at least one call of lsqfun. So, unless the residuals and their derivatives can be evaluated very quickly, the run time will be dominated by the time spent in lsqfun.
Ideally, the problem should be scaled so that, at the solution, F(x) and the corresponding values of the xj are each in the range (−1,+1), and so that at points one unit away from the solution F(x) differs from its value at the solution by approximately one unit. This will usually imply that the Hessian matrix of F(x) at the solution is well-conditioned. It is unlikely that you will be able to follow these recommendations very closely, but it is worth trying (by guesswork), as sensible scaling will reduce the difficulty of the minimization problem, so that e04gbf will take less computer time.
When the sum of squares represents the goodness-of-fit of a nonlinear model to observed data, elements of the variance-covariance matrix of the estimated regression coefficients can be computed by a subsequent call to e04ycf, using information returned in the arrays s and v. See e04ycf for further details.

10 Example

This example finds least squares estimates of x1, x2 and x3 in the model
y=x1+t1x2t2+x3t3  
using the 15 sets of data given in the following table.
y t1 t2 t3 0.14 1.0 15.0 1.0 0.18 2.0 14.0 2.0 0.22 3.0 13.0 3.0 0.25 4.0 12.0 4.0 0.29 5.0 11.0 5.0 0.32 6.0 10.0 6.0 0.35 7.0 9.0 7.0 0.39 8.0 8.0 8.0 0.37 9.0 7.0 7.0 0.58 10.0 6.0 6.0 0.73 11.0 5.0 5.0 0.96 12.0 4.0 4.0 1.34 13.0 3.0 3.0 2.10 14.0 2.0 2.0 4.39 15.0 1.0 1.0  
Before calling e04gbf, the program calls e04yaf to check lsqfun. It uses (0.5,1.0,1.5) as the initial guess at the position of the minimum.

10.1 Program Text

Program Text (e04gbfe.f90)

10.2 Program Data

Program Data (e04gbfe.d)

10.3 Program Results

Program Results (e04gbfe.r)