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_func_rcomm (c05qd)

Purpose

nag_roots_sys_func_rcomm (c05qd) is a comprehensive reverse communication function that finds a solution of a system of nonlinear equations by a modification of the Powell hybrid method.

Syntax

[irevcm, x, fvec, diag, fjac, r, qtf, iwsav, rwsav, ifail] = c05qd(irevcm, x, fvec, ml, mu, mode, diag, fjac, r, qtf, iwsav, rwsav, 'n', n, 'xtol', xtol, 'epsfcn', epsfcn, 'factor', factor)
[irevcm, x, fvec, diag, fjac, r, qtf, iwsav, rwsav, ifail] = nag_roots_sys_func_rcomm(irevcm, x, fvec, ml, mu, mode, diag, fjac, r, qtf, iwsav, rwsav, 'n', n, 'xtol', xtol, 'epsfcn', epsfcn, 'factor', factor)

Description

The system of equations is defined as:
fi (x1,x2,,xn) = 0 ,   i = 1, 2, , n .
fi (x1,x2,,xn) = 0 ,   i= 1, 2, , n .
nag_roots_sys_func_rcomm (c05qd) is based on the MINPACK routine HYBRD (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. At the starting point, the Jacobian is approximated by forward differences, but these are not used again until the rank-1 method fails to produce satisfactory progress. 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 parameter irevcm. Between intermediate exits and re-entries, all parameters other than fvec must remain unchanged.

Compulsory Input Parameters

1:     irevcm – int64int32nag_int scalar
On initial entry: must have the value 00.
Constraint: irevcm = 0irevcm=0, 11 or 22.
2:     x(n) – double array
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
On initial entry: an initial guess at the solution vector.
3:     fvec(n) – double array
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
On initial entry: need not be set.
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
On intermediate re-entry: if irevcm = 1 irevcm=1 , fvec must not be changed.
If irevcm = 2 irevcm=2 , fvec must be set to the values of the functions computed at the current point x.
4:     ml – int64int32nag_int scalar
On initial entry: the number of subdiagonals within the band of the Jacobian matrix. (If the Jacobian is not banded, or you are unsure, set ml = n1 ml=n-1 .)
Constraint: ml0 ml0 .
5:     mu – int64int32nag_int scalar
On initial entry: the number of superdiagonals within the band of the Jacobian matrix. (If the Jacobian is not banded, or you are unsure, set mu = n1 mu=n-1 .)
Constraint: mu0 mu0 .
6:     mode – int64int32nag_int scalar
On initial entry: indicates whether or not you have provided scaling factors in diag.
If mode = 2mode=2 the scaling must have been supplied in diag.
Otherwise, if mode = 1mode=1, the variables will be scaled internally.
Constraint: mode = 1mode=1 or 22.
7:     diag(n) – double array
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
If mode = 2mode=2, diag must contain multiplicative scale factors for the variables.
If mode = 1mode=1, diag need not be set.
Constraint: if mode = 2mode=2, diag(i) > 0.0 diagi>0.0 , for i = 1,2,,ni=1,2,,n.
8:     fjac(n,n) – double array
n, the first dimension of the array, must satisfy the constraint n > 0 n>0 .
On initial entry: need not be set.
9:     r(n × (n + 1) / 2n×(n+1)/2) – double array
On initial entry: need not be set.
10:   qtf(n) – double array
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
On initial entry: need not be set.
11:   iwsav(1717) – int64int32nag_int array
12:   rwsav(4 × n + 104×n+10) – double array
The arrays iwsav and rwsav must not be altered between calls to nag_roots_sys_func_rcomm (c05qd).

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.)
nn, the number of equations.
Constraint: n > 0 n>0 .
2:     xtol – double scalar
On initial entry: the accuracy in x to which the solution is required.
Suggested value: sqrt(ε)ε, where εε is the machine precision returned by nag_machine_precision (x02aj).
Default: sqrt(machine precision) machine precision
Constraint: xtol0.0 xtol0.0 .
3:     epsfcn – double scalar
On initial entry: the order of the largest relative error in the functions. It is used in determining a suitable step for a forward difference approximation to the Jacobian. If epsfcn is less than machine precision (returned by nag_machine_precision (x02aj)) then machine precision is used. Consequently a value of 0.00.0 will often be suitable.
Default: 0.00.0
4:     factor – double scalar
On initial entry: a quantity to be used in determining the initial step bound. In most cases, factor should lie between 0.10.1 and 100.0100.0. (The step bound is factor × diag × x2 factor×diag×x2  if this is nonzero; otherwise the bound is factor.)
Default: 100.0100.0
Constraint: factor > 0.0 factor>0.0 .

Input Parameters Omitted from the MATLAB Interface

None.

Output Parameters

1:     irevcm – int64int32nag_int scalar
On intermediate exit: specifies what action you must take before re-entering nag_roots_sys_func_rcomm (c05qd) with irevcm unchanged. The value of irevcm should be interpreted as follows:
irevcm = 1 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 irevcm=2
Indicates that before re-entry to nag_roots_sys_func_rcomm (c05qd), fvec must contain the function values fi(x) fi(x) .
On final exit: irevcm = 0 irevcm=0 , and the algorithm has terminated.
2:     x(n) – double array
On intermediate exit: contains the current point.
On final exit: the final estimate of the solution vector.
3:     fvec(n) – double array
On final exit: the function values at the final point, x.
4:     diag(n) – double array
The scale factors actually used (computed internally if mode = 1mode=1).
5:     fjac(n,n) – double array
On intermediate exit: must not be changed.
On final exit: the orthogonal matrix QQ produced by the QR QR  factorization of the final approximate Jacobian.
6:     r(n × (n + 1) / 2n×(n+1)/2) – double array
On intermediate exit: must not be changed.
On final exit: the upper triangular matrix RR produced by the QR QR  factorization of the final approximate Jacobian, stored row-wise.
7:     qtf(n) – double array
On intermediate exit: must not be changed.
On final exit: the vector QTf QTf .
8:     iwsav(1717) – int64int32nag_int array
9:     rwsav(4 × n + 104×n+10) – double array
10:   ifail – int64int32nag_int scalar
ifail = 0ifail=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 = 2ifail=2
Constraint: irevcm = 0irevcm=0, 11 or 22.
W ifail = 3ifail=3
No further improvement in the solution is possible.
W ifail = 4ifail=4
The iteration is not making good progress, as measured by the improvement from the last __ Jacobian evaluations.
W ifail = 5ifail=5
The iteration is not making good progress, as measured by the improvement from the last __ iterations.
  ifail = 11ifail=11
Constraint: n > 0n>0.
  ifail = 12ifail=12
Constraint: xtol0.0xtol0.0.
  ifail = 13ifail=13
Constraint: mode = 1mode=1 or 22.
  ifail = 14ifail=14
Constraint: factor > 0.0factor>0.0.
  ifail = 15ifail=15
On entry, mode = 2mode=2 and diag contained a non-positive element.
  ifail = 16ifail=16
Constraint: ml0ml0.
  ifail = 17ifail=17
Constraint: mu0mu0.
A value of ifail = 4ifail=4 or 55 may indicate that the system does not have a zero, or that the solution is very close to the origin (see Section [Accuracy]). Otherwise, rerunning nag_roots_sys_func_rcomm (c05qd) from a different starting point may avoid the region of difficulty.

Accuracy

If x^  is the true solution and DD denotes the diagonal matrix whose entries are defined by the array diag, then nag_roots_sys_func_rcomm (c05qd) tries to ensure that
D(x)2 xtol × D2 .
D (x-x^) 2 xtol × D x^ 2 .
If this condition is satisfied with xtol = 10k xtol = 10-k , then the larger components of Dx Dx  have kk significant decimal digits. There is a danger that the smaller components of Dx Dx  may have large relative errors, but the fast rate of convergence of nag_roots_sys_func_rcomm (c05qd) 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 = 3ifail=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 are reasonably well behaved. If this condition is not satisfied, then nag_roots_sys_func_rcomm (c05qd) may incorrectly indicate convergence. The validity of the answer can be checked, for example, by rerunning nag_roots_sys_func_rcomm (c05qd) with a lower value for xtol.

Further Comments

The time required by nag_roots_sys_func_rcomm (c05qd) to solve a given problem depends on nn, the behaviour of the functions, the accuracy requested and the starting point. The number of arithmetic operations executed by nag_roots_sys_func_rcomm (c05qd) to process the evaluation of functions in the main program in each exit is approximately 11.5 × n211.5×n2. The timing of nag_roots_sys_func_rcomm (c05qd) 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.
The number of function evaluations required to evaluate the Jacobian may be reduced if you can specify ml and mu accurately.

Example

function nag_roots_sys_func_rcomm_example
% The following starting values provide a rough solution.
x = -ones(9, 1);
diagnl = ones(9,1);
ml = int64(1);
mu = int64(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, diag, fjac, r, qtf, iwsav, rwsav, ifail] = ...
      nag_roots_sys_func_rcomm(irevcm, x, fvec, ml, mu, mode, diagnl, fjac, r, qtf, iwsav, rwsav);

  switch irevcm
    case {1}
      icount = icount + 1;
    case {2}
      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);
  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 {2, 3, 4}
    fprintf('\nApproximate solution\n');
    disp(x);
end
 

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


function c05qd_example
% The following starting values provide a rough solution.
x = -ones(9, 1);
diagnl = ones(9,1);
ml = int64(1);
mu = int64(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, diag, fjac, r, qtf, iwsav, rwsav, ifail] = ...
      c05qd(irevcm, x, fvec, ml, mu, mode, diagnl, fjac, r, qtf, iwsav, rwsav);

  switch irevcm
    case {1}
      icount = icount + 1;
    case {2}
      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);
  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 {2, 3, 4}
    fprintf('\nApproximate solution\n');
    disp(x);
end
 

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