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_sparsys_func_expert (c05qs)

Purpose

nag_roots_sparsys_func_expert (c05qs) is an easy-to-use function that finds a solution of a sparse system of nonlinear equations by a modification of the Powell hybrid method.

Syntax

[x, fvec, rcomm, icomm, user, ifail] = c05qs(fcn, x, init, rcomm, icomm, 'n', n, 'xtol', xtol, 'lrcomm', lrcomm, 'licomm', licomm, 'user', user)
[x, fvec, rcomm, icomm, user, ifail] = nag_roots_sparsys_func_expert(fcn, x, init, rcomm, icomm, 'n', n, 'xtol', xtol, 'lrcomm', lrcomm, 'licomm', licomm, 'user', user)

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_sparsys_func_expert (c05qs) is based on the MINPACK routine HYBRD1 (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 sparse rank-1 method of Schubert (see Schubert (1970)). At the starting point, the sparsity pattern is determined and the Jacobian is approximated by forward differences, but these are not used again until the rank-1 method fails to produce satisfactory progress. Then, the sparsity structure is used to recompute an approximation to the Jacobian by forward differences with the least number of function evaluations. The function you supply must be able to compute only the requested subset of the function values. The sparse Jacobian linear system is solved at each iteration with nag_sparse_direct_real_gen_lu (f11me) computing the Newton step. For more details see Powell (1970) and Broyden (1965).

References

Broyden C G (1965) A class of methods for solving nonlinear simultaneous equations Mathematics of Computation 19(92) 577–593
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
Schubert L K (1970) Modification of a quasi-Newton method for nonlinear equations with a sparse Jacobian Mathematics of Computation 24(109) 27–30

Parameters

Compulsory Input Parameters

1:     fcn – function handle or string containing name of m-file
fcn must return the values of the functions fifi at a point xx.
[fvec, user, iflag] = fcn(n, lindf, indf, x, user, iflag)

Input Parameters

1:     n – int64int32nag_int scalar
nn, the number of equations.
2:     lindf – int64int32nag_int scalar
lindf specifies the number of indices ii for which values of fi(x)fi(x) must be computed.
3:     indf(lindf) – int64int32nag_int array
indf specifies the indices ii for which values of fi(x)fi(x) must be computed. The indices are specified in strictly ascending order.
4:     x(n) – double array
The components of the point xx at which the functions must be evaluated. x(i)xi contains the coordinate xixi.
5:     user – Any MATLAB object
fcn is called from nag_roots_sparsys_func_expert (c05qs) with the object supplied to nag_roots_sparsys_func_expert (c05qs).
6:     iflag – int64int32nag_int scalar
iflag > 0 iflag>0 .

Output Parameters

1:     fvec(n) – double array
fvec(i)fveci must contain the function values fi(x)fi(x), for all indices ii in indf.
2:     user – Any MATLAB object
3:     iflag – int64int32nag_int scalar
In general, iflag should not be reset by fcn. If, however, you wish to terminate execution (perhaps because some illegal point x has been reached), then iflag should be set to a negative integer.
2:     x(n) – double array
n, the dimension of the array, must satisfy the constraint n > 0 n>0 .
An initial guess at the solution vector. x(i)xi must contain the coordinate xixi.
3:     init – logical scalar
init must be set to true to indicate that this is the first time nag_roots_sparsys_func_expert (c05qs) is called for this specific problem. nag_roots_sparsys_func_expert (c05qs) then computes the dense Jacobian and detects and stores its sparsity pattern (in rcomm and icomm) before proceeding with the iterations. This is noticeably time consuming when n is large. If not enough storage has been provided for rcomm or icomm, nag_roots_sparsys_func_expert (c05qs) will fail. On exit with ifail = 0ifail=0, 22, 33 or 44, icomm(1)icomm1 contains nnznnz, the number of nonzero entries found in the Jacobian. On subsequent calls, init can be set to false if the problem has a Jacobian of the same sparsity pattern. In that case, the computation time required for the detection of the sparsity pattern will be smaller.
4:     rcomm(lrcomm) – double array
lrcomm, the dimension of the array, must satisfy the constraint lrcomm12 + nnzlrcomm12+nnz where nnznnz is the number of nonzero entries in the Jacobian, as computed by nag_roots_sparsys_func_expert (c05qs).
rcomm must not be altered between successive calls to nag_roots_sparsys_func_expert (c05qs).
5:     icomm(licomm) – int64int32nag_int array
licomm, the dimension of the array, must satisfy the constraint licomm8 × n + 19 + nnzlicomm8×n+19+nnz where nnznnz is the number of nonzero entries in the Jacobian, as computed by nag_roots_sparsys_func_expert (c05qs).
If ifail = 0ifail=0, 22, 33 or 44 on exit, icomm(1)icomm1 contains nnznnz where nnznnz is the number of nonzero entries in the Jacobian.
icomm must not be altered between successive calls to nag_roots_sparsys_func_expert (c05qs).

Optional Input Parameters

1:     n – int64int32nag_int scalar
Default: The dimension of the array x.
nn, the number of equations.
Constraint: n > 0 n>0 .
2:     xtol – double scalar
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:     lrcomm – int64int32nag_int scalar
Default: The dimension of the array rcomm.
The dimension of the array rcomm as declared in the (sub)program from which nag_roots_sparsys_func_expert (c05qs) is called.
Constraint: lrcomm12 + nnzlrcomm12+nnz where nnznnz is the number of nonzero entries in the Jacobian, as computed by nag_roots_sparsys_func_expert (c05qs).
4:     licomm – int64int32nag_int scalar
Default: The dimension of the array icomm.
The dimension of the array icomm as declared in the (sub)program from which nag_roots_sparsys_func_expert (c05qs) is called.
Constraint: licomm8 × n + 19 + nnzlicomm8×n+19+nnz where nnznnz is the number of nonzero entries in the Jacobian, as computed by nag_roots_sparsys_func_expert (c05qs).
5:     user – Any MATLAB object
user is not used by nag_roots_sparsys_func_expert (c05qs), but is passed to fcn. 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

iuser ruser

Output Parameters

1:     x(n) – double array
The final estimate of the solution vector.
2:     fvec(n) – double array
The function values at the final point returned in x. fvec(i)fveci contains the function values fifi.
3:     rcomm(lrcomm) – double array
4:     icomm(licomm) – int64int32nag_int array
icomm must not be altered between successive calls to nag_roots_sparsys_func_expert (c05qs).
5:     user – Any MATLAB object
6:     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.

W ifail = 2ifail=2
There have been at least 200 × (n + 1) 200 × (n+1)  calls to fcn. Consider setting init = falseinit=false and restarting the calculation from the point held in x.
W ifail = 3ifail=3
No further improvement in the solution is possible.
W ifail = 4ifail=4
The iteration is not making good progress.
W ifail = 5ifail=5
iflag was set negative in fcn.
  ifail = 6ifail=6
lrcomm is too small.
  ifail = 7ifail=7
licomm is too small.
  ifail = 9ifail=9
An internal error has occurred. Code
  ifail = 11ifail=11
Constraint: n > 0n>0.
  ifail = 12ifail=12
Constraint: xtol0.0xtol0.0.
  ifail = 999ifail=-999
Dynamic memory allocation failed.

Accuracy

If x^  is the true solution, nag_roots_sparsys_func_expert (c05qs) tries to ensure that
x2 xtol × 2 .
x-x^ 2 xtol × x^ 2 .
If this condition is satisfied with xtol = 10k xtol = 10-k , then the larger components of xx have kk significant decimal digits. There is a danger that the smaller components of xx may have large relative errors, but the fast rate of convergence of nag_roots_sparsys_func_expert (c05qs) 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_sparsys_func_expert (c05qs) may incorrectly indicate convergence. The validity of the answer can be checked, for example, by rerunning nag_roots_sparsys_func_expert (c05qs) with a lower value for xtol.

Further Comments

Local workspace arrays of fixed lengths are allocated internally by nag_roots_sparsys_func_expert (c05qs). The total size of these arrays amounts to 8 × n + 2 × q8×n+2×q double elements and 10 × n + 2 × q + 510×n+2×q+5 integer elements where the integer qq is bounded by 8 × nnz8×nnz and n2n2 and depends on the sparsity pattern of the Jacobian.
The time required by nag_roots_sparsys_func_expert (c05qs) 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_sparsys_func_expert (c05qs) to process each evaluation of the functions depends on the number of nonzero entries in the Jacobian. The timing of nag_roots_sparsys_func_expert (c05qs) is strongly influenced by the time spent evaluating the functions.
When init is true, the dense Jacobian is first evaluated and that will take time proportional to n2n2.
Ideally the problem should be scaled so that, at the solution, the function values are of comparable magnitude.

Example

function nag_roots_sparsys_func_expert_example
% The following starting values provide a rough solution.
x = -ones(9, 1);
rcomm = zeros(39, 1);
icomm = zeros(118, 1, 'int64');
for i=0:1

  user = i; % Perturb the system?
  init = (i==0);

  [x, fvec, rcomm, icomm, user, ifail] = ...
     nag_roots_sparsys_func_expert(@fcn, x, init, rcomm, icomm,'user', user);

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

end



function [fvec, user, iflag] = fcn(n, lindf, indf, x, user, iflag)
  fvec = zeros(n, 1);
  ind = 1;
  iflag = int64(0);
  alpha = (1/2)^7;
  theta = user*alpha;
  for i = 1:double(n)
    if indf(ind) ~= i
      continue;
    end
    fvec(i) = (3-(2+theta)*x(i))*x(i) + 1;
    if (i > 1)
      fvec(i) = fvec(i) - x(i-1);
    end
    if (i < n)
      fvec(i) = fvec(i) - 2*x(i+1);
    end
    ind = ind + 1;
    if (ind > lindf)
     break;
    end
  end
 

Final 2-norm of the residuals =   1.7592e-09

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


Final 2-norm of the residuals =   2.6329e-13

Final approximate solution
   -0.5697
   -0.6804
   -0.7004
   -0.7029
   -0.7000
   -0.6906
   -0.6646
   -0.5951
   -0.4159


function c05qs_example
% The following starting values provide a rough solution.
x = -ones(9, 1);
rcomm = zeros(39, 1);
icomm = zeros(118, 1, 'int64');
for i=0:1

  user = i; % Perturb the system?
  init = (i==0);

  [x, fvec, rcomm, icomm, user, ifail] = ...
     c05qs(@fcn, x, init, rcomm, icomm,'user', user);

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

end



function [fvec, user, iflag] = fcn(n, lindf, indf, x, user, iflag)
  fvec = zeros(n, 1);
  ind = 1;
  iflag = int64(0);
  alpha = (1/2)^7;
  theta = user*alpha;
  for i = 1:double(n)
    if indf(ind) ~= i
      continue;
    end
    fvec(i) = (3-(2+theta)*x(i))*x(i) + 1;
    if (i > 1)
      fvec(i) = fvec(i) - x(i-1);
    end
    if (i < n)
      fvec(i) = fvec(i) - 2*x(i+1);
    end
    ind = ind + 1;
    if (ind > lindf)
     break;
    end
  end
 

Final 2-norm of the residuals =   1.7592e-09

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


Final 2-norm of the residuals =   2.6329e-13

Final approximate solution
   -0.5697
   -0.6804
   -0.7004
   -0.7029
   -0.7000
   -0.6906
   -0.6646
   -0.5951
   -0.4159



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