hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_opt_lsq_uncon_quasi_deriv_comp (e04gb)

Purpose

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

Syntax

[x, fsumsq, fvec, fjac, s, v, niter, nf, user, ifail] = e04gb(m, lsqlin, lsqfun, lsqmon, maxcal, eta, xtol, x, 'n', n, 'iprint', iprint, 'stepmx', stepmx, 'user', user)
[x, fsumsq, fvec, fjac, s, v, niter, nf, user, ifail] = nag_opt_lsq_uncon_quasi_deriv_comp(m, lsqlin, lsqfun, lsqmon, maxcal, eta, xtol, x, 'n', n, 'iprint', iprint, 'stepmx', stepmx, 'user', user)
Note: the interface to this routine has changed since earlier releases of the toolbox:
Mark 22: liw, lw have been removed from the interface
Mark 24: drop w and iw, introduce user in main routine and all call-backs
.

Description

nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) is essentially identical to the function LSQFDQ in the NPL Algorithms Library. It is applicable to problems of the form:
m
MinimizeF(x) = [fi(x)]2
i = 1
MinimizeF(x)=i=1m[fi(x)]2
where x = (x1,x2,,xn)T x = (x1,x2,,xn)T  and mnmn. (The functions fi(x)fi(x) are often referred to as ‘residuals’.)
You must supply a function to calculate the values of the fi(x)fi(x) and their first derivatives (fi)/(xj) fi xj  at any point xx.
From a starting point x(1)x (1)  supplied by you, the function generates a sequence of points x(2),x(3),x (2) ,x (3) ,, which is intended to converge to a local minimum of F(x)F(x). The sequence of points is given by
x(k + 1) = x(k) + α(k)p(k)
x (k+1) =x (k) +α(k)p (k)
where the vector p(k)p (k)  is a direction of search, and α(k)α(k) is chosen such that F(x(k) + α(k)p(k))F(x (k) +α(k)p (k) ) is approximately a minimum with respect to α(k)α(k).
The vector p(k)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)p (k)  is the Gauss–Newton direction; otherwise the second derivatives of the fi(x)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.

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

Parameters

Compulsory Input Parameters

1:     m – int64int32nag_int scalar
The number mm of residuals, fi(x)fi(x), and the number nn of variables, xjxj.
Constraint: 1nm1nm.
2:     lsqlin – function handle or string containing name of m-file
lsqlin enables you to specify whether the linear minimizations (i.e., minimizations of F(x(k) + α(k)p(k))F(x (k) +α (k) p (k) ) with respect to α(k)α(k)) are to be performed by a function which just requires the evaluation of the fi(x)fi(x) (nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_fun (e04fcv)), or by a function which also requires the first derivatives of the fi(x)fi(x) (nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_deriv (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 nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) should be called with string 'e04hev'>. However, if the evaluation of the derivatives takes more than about 44 times as long as the evaluation of the residuals, then nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_fun (e04fcv) will usually be preferable. If in doubt, use as it is slightly more robust.
3:     lsqfun – function handle or string containing name of m-file
lsqfun must calculate the vector of values fi(x)fi(x) and Jacobian matrix of first derivatives (fi)/(xj) fi xj  at any point xx. (However, if you do not wish to calculate the residuals or first derivatives at a particular xx, there is the option of setting a parameter to cause nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) to terminate immediately.)
[iflag, fvec, fjac, user] = lsqfun(iflag, m, n, xc, ldfjac, user)

Input Parameters

1:     iflag – int64int32nag_int scalar
Will be set to 00, 11 or 22.
iflag = 0iflag=0
Indicates that only the residuals need to be evaluated
iflag = 1iflag=1
Indicates that only the Jacobian matrix needs to be evaluated
iflag = 2iflag=2
Indicates that both the residuals and the Jacobian matrix must be calculated.
If nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_deriv (e04hev) is used as nag_opt_lsq_uncon_quasi_deriv_comp (e04gb)'s lsqlin, lsqfun will always be called with iflag set to 22.
2:     m – int64int32nag_int scalar
mm, the number of residuals.
3:     n – int64int32nag_int scalar
nn, the number of variables.
4:     xc(n) – double array
The point xx at which the values of the fifi and the (fi)/(xj) fi xj  are required.
5:     ldfjac – int64int32nag_int scalar
The first dimension of the array fjac as declared in the (sub)program from which nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) is called.
6:     user – Any MATLAB object
lsqfun is called from nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) with the object supplied to nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).

Output Parameters

1:     iflag – int64int32nag_int scalar
If it is not possible to evaluate the fi(x)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 nag_opt_lsq_uncon_quasi_deriv_comp (e04gb). nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) will then terminate immediately, with ifail set to your setting of iflag.
2:     fvec(m) – double array
Unless iflag = 1iflag=1 on entry, or iflag is reset to a negative number, then fvec(i)fveci must contain the value of fifi at the point xx, for i = 1,2,,mi=1,2,,m.
3:     fjac(ldfjac,n) – double array
ldfjacmldfjacm.
Unless iflag = 0iflag=0 on entry, or iflag is reset to a negative number, then fjac(i,j)fjacij must contain the value of (fi)/(xj) fi xj at the point xx, for i = 1,2,,mi=1,2,,m and j = 1,2,,nj=1,2,,n.
4:     user – Any MATLAB object
Note:  lsqfun should be tested separately before being used in conjunction with nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).
4:     lsqmon – function handle or string containing name of m-file
If iprint0iprint0, you must supply lsqmon which is suitable for monitoring the minimization process. lsqmon must not change the values of any of its parameters.
If iprint < 0iprint<0, the NAG Toolbox string 'e04fdz' can be used as lsqmon.
[user] = lsqmon(m, n, xc, fvec, fjac, ldfjac, s, igrade, niter, nf, user)

Input Parameters

1:     m – int64int32nag_int scalar
mm, the numbers of residuals.
2:     n – int64int32nag_int scalar
nn, the numbers of variables.
3:     xc(n) – double array
The coordinates of the current point xx.
4:     fvec(m) – double array
The values of the residuals fifi at the current point xx.
5:     fjac(ldfjac,n) – double array
fjac(i,j)fjacij contains the value of (fi)/(xj) fi xj  at the current point xx, for i = 1,2,,mi=1,2,,m and j = 1,2,,nj=1,2,,n.
6:     ldfjac – int64int32nag_int scalar
The first dimension of the array fjac as declared in the (sub)program from which nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) is called.
7:     s(n) – double array
The singular values of the current Jacobian matrix. Thus s may be useful as information about the structure of your problem.
8:     igrade – int64int32nag_int scalar
nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) 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 – int64int32nag_int scalar
The number of iterations which have been performed in nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).
10:   nf – int64int32nag_int scalar
The number of evaluations of the residuals. (If nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_deriv (e04hev) is used as lsqlin, nf is also the number of evaluations of the Jacobian matrix.)
11:   user – Any MATLAB object
lsqmon is called from nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) with the object supplied to nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).

Output Parameters

1:     user – Any MATLAB object
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)F(x) mentioned in Section [Accuracy]. It is usually helpful to also print xc, the gradient of the sum of squares, niter and nf.
5:     maxcal – int64int32nag_int scalar
Enables you to limit the number of times that lsqfun is called by nag_opt_lsq_uncon_quasi_deriv_comp (e04gb). There will be an error exit (see Section [Error Indicators and Warnings]) after maxcal calls of lsqfun.
Constraint: maxcal1maxcal1.
6:     eta – double scalar
Every iteration of nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) involves a linear minimization (i.e., minimization of F(x(k) + α(k)p(k))F(x (k) +α (k) p (k) ) with respect to α(k)α (k) ). eta specifies how accurately these linear minimizations are to be performed. The minimum with respect to α(k)α (k)  will be located more accurately for small values of eta (say, 0.010.01) than for large values (say, 0.90.9).
Although accurate linear minimizations will generally reduce the number of iterations performed by nag_opt_lsq_uncon_quasi_deriv_comp (e04gb), 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.
Constraint: 0.0eta < 1.00.0eta<1.0.
7:     xtol – double scalar
The accuracy in xx to which the solution is required.
If xtruextrue is the true value of xx at the minimum, then xsolxsol, the estimated position before a normal exit, is such that
xsolxtrue < xtol × (1.0 + xtrue),
xsol-xtrue<xtol×(1.0+xtrue),
where y = sqrt(j = 1nyj2)y=j=1nyj2. For example, if the elements of xsolxsol are not much larger than 1.01.0 in modulus and if xtol = 1.0e−5xtol=1.0e−5, then xsolxsol is usually accurate to about five decimal places. (For further details see Section [Accuracy].)
If F(x)F(x) and the variables are scaled roughly as described in Section [Further Comments] and εε is the machine precision, then a setting of order xtol = sqrt(ε)xtol=ε will usually be appropriate. If xtol is set to 0.00.0 or some positive value less than 10ε10ε, nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) will use 10ε10ε instead of xtol, since 10ε10ε is probably the smallest reasonable setting.
Constraint: xtol0.0xtol0.0.
8:     x(n) – double array
n, the dimension of the array, must satisfy the constraint 1nm1nm.
x(j)xj must be set to a guess at the jjth component of the position of the minimum, for j = 1,2,,nj=1,2,,n.

Optional Input Parameters

1:     n – int64int32nag_int scalar
Default: For n, the dimension of the array x.
The number mm of residuals, fi(x)fi(x), and the number nn of variables, xjxj.
Constraint: 1nm1nm.
2:     iprint – int64int32nag_int scalar
The frequency with which lsqmon is to be called.
iprint > 0iprint>0
lsqmon is called once every iprint iterations and just before exit from nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).
iprint = 0iprint=0
lsqmon is just called at the final point.
iprint < 0iprint<0
lsqmon is not called at all.
iprint should normally be set to a small positive number.
Default: 11
3:     stepmx – double scalar
An estimate of the Euclidean distance between the solution and the starting point supplied by you. (For maximum efficiency, a slight overestimate is preferable.)
nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) will ensure that, for each iteration,
n
(xj(k)xj(k1))2(stepmx)2
j = 1
j=1n(xj (k) -xj (k-1) )2 (stepmx) 2
where kk is the iteration number. Thus, if the problem has more than one solution, nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) 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)x (k)  entering a region where the problem is ill-behaved and can help avoid overflow in the evaluation of F(x)F(x). However, an underestimate of stepmx can lead to inefficiency.
Default: 100000.0100000.0
Constraint: stepmxxtolstepmxxtol.
4:     user – Any MATLAB object
user is not used by nag_opt_lsq_uncon_quasi_deriv_comp (e04gb), but is passed to lsqfun and lsqmon. Note that for large objects it may be more efficient to use a global variable which is accessible from the m-files than to use user.

Input Parameters Omitted from the MATLAB Interface

ldfjac ldv iw liw w lw

Output Parameters

1:     x(n) – double array
The final point x(k)x (k) . Thus, if ifail = 0ifail=0 on exit, x(j)xj is the jjth component of the estimated position of the minimum.
2:     fsumsq – double scalar
The value of F(x)F(x), the sum of squares of the residuals fi(x)fi(x), at the final point given in x.
3:     fvec(m) – double array
The value of the residual fi(x)fi(x) at the final point given in x, for i = 1,2,,mi=1,2,,m.
4:     fjac(ldfjac,n) – double array
ldfjacmldfjacm.
The value of the first derivative (fi)/(xj) fi xj evaluated at the final point given in x, for i = 1,2,,mi=1,2,,m and j = 1,2,,nj=1,2,,n.
5:     s(n) – double array
The singular values of the Jacobian matrix at the final point. Thus s may be useful as information about the structure of your problem.
6:     v(ldv,n) – double array
ldvnldvn.
The matrix VV associated with the singular value decomposition
J = USVT
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 JTJJTJ.
7:     niter – int64int32nag_int scalar
The number of iterations which have been performed in nag_opt_lsq_uncon_quasi_deriv_comp (e04gb).
8:     nf – int64int32nag_int scalar
The number of times that the residuals have been evaluated (i.e., the number of calls of lsqfun). If nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_deriv (e04hev) is used as lsqlin, nf is also the number of times that the Jacobian matrix has been evaluated.
9:     user – Any MATLAB object
10:   ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Note: nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) may return useful information for one or more of the following detected errors or warnings.
Errors or warnings detected by the function:

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

W ifail < 0ifail<0
A negative value of ifail indicates an exit from nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) because you have set iflag negative in lsqfun. The value of ifail will be the same as your setting of iflag.
  ifail = 1ifail=1
On entry,n < 1n<1,
orm < nm<n,
ormaxcal < 1maxcal<1,
oreta < 0.0eta<0.0,
oreta1.0eta1.0,
orxtol < 0.0xtol<0.0,
orstepmx < xtolstepmx<xtol,
orldfjac < mldfjac<m,
orldv < nldv<n,
orliw < 1liw<1,
orlw < 7 × n + m × n + 2 × m + n × nlw<7×n+m×n+2×m+n×n when n > 1n>1,
orlw < 9 + 3 × mlw<9+3×m when n = 1n=1.
When this exit occurs, no values will have been assigned to fsumsq, or to the elements of fvec, fjac, s or v.
  ifail = 2ifail=2
There have been maxcal calls of lsqfun. If steady reductions in the sum of squares, F(x)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)F(x) has no minimum.
W ifail = 3ifail=3
The conditions for a minimum have not all been satisfied, but a lower point could not be found. This could be because xtol has been set so small that rounding errors in the evaluation of the residuals and derivatives make attainment of the convergence conditions impossible.
  ifail = 4ifail=4
The method for computing the singular value decomposition of the Jacobian matrix has failed to converge in a reasonable number of sub-iterations. It may be worth applying nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) again starting with an initial approximation which is not too close to the point at which the failure occurred.
The values ifail = 2ifail=2, 33 or 44 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.

Accuracy

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

Further Comments

The number of iterations required depends on the number of variables, the number of residuals, the behaviour of F(x)F(x), the accuracy demanded and the distance of the starting point from the solution. The number of multiplications performed per iteration of nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) varies, but for mnmn is approximately n × m2 + O(n3)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)F(x) and the corresponding values of the xjxj are each in the range (1, + 1)(-1,+1), and so that at points one unit away from the solution F(x)F(x) differs from its value at the solution by approximately one unit. This will usually imply that the Hessian matrix of F(x)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 nag_opt_lsq_uncon_quasi_deriv_comp (e04gb) 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 nag_opt_lsq_uncon_covariance (e04yc), using information returned in the arrays s and v. See nag_opt_lsq_uncon_covariance (e04yc) for further details.

Example

function nag_opt_lsq_uncon_quasi_deriv_comp_example
m = int64(15);
maxcal = int64(150);
eta = 0.9;
xtol = 1.05418557512311e-07;
x = [0.5;
     1;
     1.5];
y=[0.14,0.18,0.22,0.25,0.29,0.32,0.35,0.39,0.37,0.58,0.73,0.96,1.34,2.10,4.39];
t = [1.0, 15.0, 1.0;
     2.0, 14.0, 2.0;
     3.0, 13.0, 3.0;
     4.0, 12.0, 4.0;
     5.0, 11.0, 5.0;
     6.0, 10.0, 6.0;
     7.0, 9.0, 7.0;
     8.0, 8.0, 8.0;
     9.0, 7.0, 7.0;
     10.0, 6.0, 6.0;
     11.0, 5.0, 5.0;
     12.0, 4.0, 4.0;
     13.0, 3.0, 3.0;
     14.0, 2.0, 2.0;
     15.0, 1.0, 1.0];
user = {y; t};
[xOut, fsumsq, fvec, fjac, s, v, niter, nf, user, ifail] = ...
  nag_opt_lsq_uncon_quasi_deriv_comp(m, ...
                         'nag_opt_lsq_uncon_quasi_deriv_comp_lsqlin_deriv', ...
                         @lsqfun, @lsqmon, maxcal, eta, xtol, x, 'user', user)

function [iflag, fvecc, fjacc, user] = lsqfun(iflag, m, n, xc, ljc, user)
  y = user{1};
  t = user{2};

  fvecc = zeros(m, 1);
  fjacc = zeros(ljc, n);

  for i = 1:double(m)
    denom = xc(2)*t(i,2) + xc(3)*t(i,3);
    fvecc(i) = xc(1) + t(i,1)/denom - y(i);
    if (iflag ~= 0)
      fjacc(i,1) = 1;
      dummy = -1/(denom*denom);
      fjacc(i,2) = t(i,1)*t(i,2)*dummy;
      fjacc(i,3) = t(i,1)*t(i,3)*dummy;
    end
  end

function [user] = lsqmon(m, n, xc, fvecc, fjacc, ljc, s, igrade, niter, nf, user)
 

xOut =

    0.0824
    1.1330
    2.3437


fsumsq =

    0.0082


fvec =

   -0.0059
   -0.0003
    0.0003
    0.0065
   -0.0008
   -0.0013
   -0.0045
   -0.0200
    0.0822
   -0.0182
   -0.0148
   -0.0147
   -0.0112
   -0.0042
    0.0068


fjac =

    1.0000   -0.0401   -0.0027
    1.0000   -0.0663   -0.0095
    1.0000   -0.0824   -0.0190
    1.0000   -0.0910   -0.0303
    1.0000   -0.0941   -0.0428
    1.0000   -0.0931   -0.0558
    1.0000   -0.0890   -0.0692
    1.0000   -0.0827   -0.0827
    1.0000   -0.1064   -0.1064
    1.0000   -0.1379   -0.1379
    1.0000   -0.1820   -0.1820
    1.0000   -0.2482   -0.2482
    1.0000   -0.3585   -0.3585
    1.0000   -0.5791   -0.5791
    1.0000   -1.2409   -1.2409


s =

    4.0965
    1.5950
    0.0613


v =

   -0.9354    0.3530    0.0214
    0.2592    0.6432    0.7205
    0.2405    0.6795   -0.6932


niter =

                    6


nf =

                    7


user = 

    [ 1x15 double]
    [15x3  double]


ifail =

                    0


function e04gb_example
m = int64(15);
maxcal = int64(150);
eta = 0.9;
xtol = 1.05418557512311e-07;
x = [0.5;
     1;
     1.5];
y=[0.14,0.18,0.22,0.25,0.29,0.32,0.35,0.39,0.37,0.58,0.73,0.96,1.34,2.10,4.39];
t = [1.0, 15.0, 1.0;
     2.0, 14.0, 2.0;
     3.0, 13.0, 3.0;
     4.0, 12.0, 4.0;
     5.0, 11.0, 5.0;
     6.0, 10.0, 6.0;
     7.0, 9.0, 7.0;
     8.0, 8.0, 8.0;
     9.0, 7.0, 7.0;
     10.0, 6.0, 6.0;
     11.0, 5.0, 5.0;
     12.0, 4.0, 4.0;
     13.0, 3.0, 3.0;
     14.0, 2.0, 2.0;
     15.0, 1.0, 1.0];
user = {y; t};
[xOut, fsumsq, fvec, fjac, s, v, niter, nf, user, ifail] = ...
     e04gb(m, 'e04hev', @lsqfun, @lsqmon, maxcal, eta, xtol, x, 'user', user)

function [iflag, fvecc, fjacc, user] = lsqfun(iflag, m, n, xc, ljc, user)
  y = user{1};
  t = user{2};

  fvecc = zeros(m, 1);
  fjacc = zeros(ljc, n);

  for i = 1:double(m)
    denom = xc(2)*t(i,2) + xc(3)*t(i,3);
    fvecc(i) = xc(1) + t(i,1)/denom - y(i);
    if (iflag ~= 0)
      fjacc(i,1) = 1;
      dummy = -1/(denom*denom);
      fjacc(i,2) = t(i,1)*t(i,2)*dummy;
      fjacc(i,3) = t(i,1)*t(i,3)*dummy;
    end
  end

function [user] = lsqmon(m, n, xc, fvecc, fjacc, ljc, s, igrade, niter, nf, user)
 

xOut =

    0.0824
    1.1330
    2.3437


fsumsq =

    0.0082


fvec =

   -0.0059
   -0.0003
    0.0003
    0.0065
   -0.0008
   -0.0013
   -0.0045
   -0.0200
    0.0822
   -0.0182
   -0.0148
   -0.0147
   -0.0112
   -0.0042
    0.0068


fjac =

    1.0000   -0.0401   -0.0027
    1.0000   -0.0663   -0.0095
    1.0000   -0.0824   -0.0190
    1.0000   -0.0910   -0.0303
    1.0000   -0.0941   -0.0428
    1.0000   -0.0931   -0.0558
    1.0000   -0.0890   -0.0692
    1.0000   -0.0827   -0.0827
    1.0000   -0.1064   -0.1064
    1.0000   -0.1379   -0.1379
    1.0000   -0.1820   -0.1820
    1.0000   -0.2482   -0.2482
    1.0000   -0.3585   -0.3585
    1.0000   -0.5791   -0.5791
    1.0000   -1.2409   -1.2409


s =

    4.0965
    1.5950
    0.0613


v =

   -0.9354    0.3530    0.0214
    0.2592    0.6432    0.7205
    0.2405    0.6795   -0.6932


niter =

                    6


nf =

                    7


user = 

    [ 1x15 double]
    [15x3  double]


ifail =

                    0



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013