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_roots_sys_deriv_rcomm (c05rd)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_roots_sys_deriv_rcomm (c05rd) is a comprehensive reverse communication function that finds a solution of a system of nonlinear equations by a modification of the Powell hybrid method. You must provide the Jacobian.

Syntax

[irevcm, x, fvec, fjac, diag, r, qtf, iwsav, rwsav, ifail] = c05rd(irevcm, x, fvec, fjac, mode, diag, r, qtf, iwsav, rwsav, 'n', n, 'xtol', xtol, 'factor', factor)
[irevcm, x, fvec, fjac, diag, r, qtf, iwsav, rwsav, ifail] = nag_roots_sys_deriv_rcomm(irevcm, x, fvec, fjac, mode, diag, r, qtf, iwsav, rwsav, 'n', n, 'xtol', xtol, 'factor', factor)

Description

The system of equations is defined as:
fi x1,x2,,xn = 0 ,   i= 1, 2, , n .  
nag_roots_sys_deriv_rcomm (c05rd) is based on the MINPACK routine HYBRJ (see Moré et al. (1980)). It chooses the correction at each step as a convex combination of the Newton and scaled gradient directions. The Jacobian is updated by the rank-1 method of Broyden. For more details see Powell (1970).

References

Moré J J, Garbow B S and Hillstrom K E (1980) User guide for MINPACK-1 Technical Report ANL-80-74 Argonne National Laboratory
Powell M J D (1970) A hybrid method for nonlinear algebraic equations Numerical Methods for Nonlinear Algebraic Equations (ed P Rabinowitz) Gordon and Breach

Parameters

Note: this function uses reverse communication. Its use involves an initial entry, intermediate exits and re-entries, and a final exit, as indicated by the argument irevcm. Between intermediate exits and re-entries, all arguments other than fvec and fjac must remain unchanged.

Compulsory Input Parameters

1:     irevcm int64int32nag_int scalar
On initial entry: must have the value 0.
Constraint: irevcm=0, 1, 2 or 3.
2:     xn – double array
On initial entry: an initial guess at the solution vector.
3:     fvecn – double array
On initial entry: need not be set.
On intermediate re-entry: if irevcm2 , fvec must not be changed.
If irevcm=2 , fvec must be set to the values of the functions computed at the current point x.
4:     fjacnn – double array
On initial entry: need not be set.
On intermediate re-entry: if irevcm3 , fjac must not be changed.
If irevcm=3 , fjacij  must contain the value of fi xj  at the point x, for i=1,2,,n and j=1,2,,n.
5:     mode int64int32nag_int scalar
On initial entry: indicates whether or not you have provided scaling factors in diag.
If mode=2 the scaling must have been supplied in diag.
Otherwise, if mode=1, the variables will be scaled internally.
Constraint: mode=1 or 2.
6:     diagn – double array
On initial entry: if mode=2, diag must contain multiplicative scale factors for the variables.
If mode=1, diag need not be set.
Constraint: if mode=2, diagi>0.0 , for i=1,2,,n.
7:     rn×n+1/2 – double array
On initial entry: need not be set.
8:     qtfn – double array
On initial entry: need not be set.
9:     iwsav17 int64int32nag_int array
10:   rwsav4×n+10 – double array
The arrays iwsav and rwsav must not be altered between calls to nag_roots_sys_deriv_rcomm (c05rd).

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays x, fvec, diag, qtf and the first dimension of the array fjac and the second dimension of the array fjac. (An error is raised if these dimensions are not equal.)
n, the number of equations.
Constraint: n>0 .
2:     xtol – double scalar
Suggested value: ε, where ε is the machine precision returned by nag_machine_precision (x02aj).
Default: machine precision
On initial entry: the accuracy in x to which the solution is required.
Constraint: xtol0.0 .
3:     factor – double scalar
Default: 100.0
On initial entry: a quantity to be used in determining the initial step bound. In most cases, factor should lie between 0.1 and 100.0. (The step bound is factor×diag×x2  if this is nonzero; otherwise the bound is factor.)
Constraint: factor>0.0 .

Output Parameters

1:     irevcm int64int32nag_int scalar
On intermediate exit: specifies what action you must take before re-entering nag_roots_sys_deriv_rcomm (c05rd) with irevcm unchanged. The value of irevcm should be interpreted as follows:
irevcm=1
Indicates the start of a new iteration. No action is required by you, but x and fvec are available for printing.
irevcm=2
Indicates that before re-entry to nag_roots_sys_deriv_rcomm (c05rd), fvec must contain the function values fix .
irevcm=3
Indicates that before re-entry to nag_roots_sys_deriv_rcomm (c05rd), fjacij  must contain the value of fi xj  at the point x, for i=1,2,,n and j=1,2,,n.
On final exit: irevcm=0 , and the algorithm has terminated.
2:     xn – double array
On intermediate exit: contains the current point.
On final exit: the final estimate of the solution vector.
3:     fvecn – double array
On final exit: the function values at the final point, x.
4:     fjacnn – double array
On final exit: the orthogonal matrix Q produced by the QR  factorization of the final approximate Jacobian.
5:     diagn – double array
On intermediate exit: diag must not be changed.
On final exit: the scale factors actually used (computed internally if mode=1).
6:     rn×n+1/2 – double array
On intermediate exit: must not be changed.
On final exit: the upper triangular matrix R produced by the QR  factorization of the final approximate Jacobian, stored row-wise.
7:     qtfn – double array
On intermediate exit: must not be changed.
On final exit: the vector QTf .
8:     iwsav17 int64int32nag_int array
9:     rwsav4×n+10 – double array
10:   ifail int64int32nag_int scalar
On final exit: ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and 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.

   ifail=2
Constraint: irevcm=0, 1, 2 or 3.
W  ifail=3
No further improvement in the solution is possible.
W  ifail=4
The iteration is not making good progress, as measured by the improvement from the last _ Jacobian evaluations. This failure exit may indicate that the system does not have a zero, or that the solution is very close to the origin (see Accuracy). Otherwise, rerunning nag_roots_sys_deriv_rcomm (c05rd) from a different starting point may avoid the region of difficulty.
W  ifail=5
The iteration is not making good progress, as measured by the improvement from the last _ iterations. This failure exit may indicate that the system does not have a zero, or that the solution is very close to the origin (see Accuracy). Otherwise, rerunning nag_roots_sys_deriv_rcomm (c05rd) from a different starting point may avoid the region of difficulty.
   ifail=11
Constraint: n>0.
   ifail=12
Constraint: xtol0.0.
   ifail=13
Constraint: mode=1 or 2.
   ifail=14
Constraint: factor>0.0.
   ifail=15
On entry, mode=2 and diag contained a non-positive element.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

If x^  is the true solution and D denotes the diagonal matrix whose entries are defined by the array diag, then nag_roots_sys_deriv_rcomm (c05rd) tries to ensure that
D x-x^ 2 xtol × D x^ 2 .  
If this condition is satisfied with xtol = 10-k , then the larger components of Dx  have k significant decimal digits. There is a danger that the smaller components of Dx  may have large relative errors, but the fast rate of convergence of nag_roots_sys_deriv_rcomm (c05rd) usually obviates this possibility.
If xtol is less than machine precision and the above test is satisfied with the machine precision in place of xtol, then the function exits with ifail=3.
Note:  this convergence test is based purely on relative error, and may not indicate convergence if the solution is very close to the origin.
The convergence test assumes that the functions and the Jacobian are coded consistently and that the functions are reasonably well behaved. If these conditions are not satisfied, then nag_roots_sys_deriv_rcomm (c05rd) may incorrectly indicate convergence. The coding of the Jacobian can be checked using nag_roots_sys_deriv_check (c05zd). If the Jacobian is coded correctly, then the validity of the answer can be checked by rerunning nag_roots_sys_deriv_rcomm (c05rd) with a lower value for xtol.

Further Comments

The time required by nag_roots_sys_deriv_rcomm (c05rd) to solve a given problem depends on n, the behaviour of the functions, the accuracy requested and the starting point. The number of arithmetic operations executed by nag_roots_sys_deriv_rcomm (c05rd) is approximately 11.5×n2  to process each evaluation of the functions and approximately 1.3×n3  to process each evaluation of the Jacobian. The timing of nag_roots_sys_deriv_rcomm (c05rd) is strongly influenced by the time spent evaluating the functions.
Ideally the problem should be scaled so that, at the solution, the function values are of comparable magnitude.

Example

This example determines the values x1 , , x9  which satisfy the tridiagonal equations:
3-2x1x1-2x2 = -1, -xi-1+3-2xixi-2xi+1 = -1,  i=2,3,,8 -x8+3-2x9x9 = -1.  
function c05rd_example


fprintf('c05rd example results\n\n');

% The following starting values provide a rough solution.
x = -ones(9, 1);
diagnl = ones(9,1);
mode = int64(2);
irevcm = int64(0);
fjac = zeros(9, 9);
fvec = zeros(9, 1);
qtf = zeros(9, 1);
r = zeros(45, 1);
rwsav = zeros(46, 1);
iwsav = zeros(17, 1, 'int64');

icount = int64(0);
first = true;
while (irevcm ~= 0 || first )
  first = false;
  [irevcm, x, fvec, fjac, diag, r, qtf, iwsav, rwsav, ifail] = ...
      c05rd(irevcm, x, fvec, fjac, mode, diagnl, r, qtf, iwsav, rwsav);

  switch irevcm
    case {1}
      icount = icount + 1;
    case {2}
      % Evaluate functions at given point
      fvec(1:9) = (3.0-2.0.*x).*x + 1.0;
      fvec(2:9) = fvec(2:9) - x(1:8);
      fvec(1:8) = fvec(1:8) - 2.0.*x(2:9);
    case {3}
      % Evaluate Jacobian at current point
      fjac = zeros(9, 9);
      for k = 1:9
        fjac(k,k) = 3 - 4*x(k);
        if (k ~= 1)
          fjac(k,k-1) = -1;
        end
        if (k ~= 9)
          fjac(k,k+1) = -2;
        end
      end
  end
end

switch ifail
  case {0}
    fprintf('\nFinal 2-norm of the residuals after %d iterations = %12.4e\n',...
      icount, norm(fvec));
    fprintf('\nFinal approximate solution\n');
    disp(x);
  case {3, 4, 5}
    fprintf('\nApproximate solution\n');
    disp(x);
end


c05rd example results


Final 2-norm of the residuals after 11 iterations =   1.1926e-08

Final approximate solution
   -0.5707
   -0.6816
   -0.7017
   -0.7042
   -0.7014
   -0.6919
   -0.6658
   -0.5960
   -0.4164


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–2015