| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| integer, | intent(in) | :: | n | Dimension of the problem |
||
| real(kind=r8), | intent(inout), | dimension(n) | :: | x | Vector of length n, initial guess on input and final approximate solution on output |
|
| procedure(f_intr) | :: | f | User-supplied subroutine for evaluating the function the zero of which is sought. |
|||
| procedure(jacv_intr) | :: | jacv | User-supplied subroutine for optionally evaluating or , where is the Jacobian of and is a right preconditioning operator. If neither analytic evaluations nor right preconditioning is used, this can be a dummy subroutine; if right preconditioning is used but not analytic evaluations, this need only evaluate . |
|||
| real(kind=r8), | intent(in) | :: | ftol | Stopping tolerance of the f-norm |
||
| real(kind=r8), | intent(in) | :: | stptol | Stopping tolerance of the step-length |
||
| integer, | intent(in), | dimension(10) | :: | input | Array containing various user-specified inputs; see above |
|
| integer, | intent(out), | dimension(6) | :: | info | Array containing various outputs; see above |
|
| real(kind=r8), | intent(inout), | dimension(*) | :: | rwork | Work array with length depending on the solver used as follows:
|
|
| real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | Parameter/work array passed to the |
|
| integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array passed to the |
|
| integer, | intent(out) | :: | iterm | Termination flag. Values have the following meanings:
|
||
| procedure(dinpr_intr) | :: | dinpr | User-supplied function for calculating vector inner products.
This has the same interace as the BLAS routine
ddot.
If the Euclidean inner product is desired then user can link
to a local BLAS library and provide the name |
|||
| procedure(dnorm_intr) | :: | dnorm | User-supplied function for calculating vector norms. This has the same interface as the BLAS routine dnrm2; if the Euclidean norm is desired the user can link to a local BLAS library and provide the name dnrm2 to nitsol. dnorm must be declared as an external function that returns a double precision value in the calling program. |
An explicit interface to the nitsol Newton iterative nonlinear solver.
The input array argument allows the user to specify various
options. It should be declared an integer vector of length 11
[sic.] in the calling program. To specify an option, set the
appropriate input component to the desired value according to
the specifications below.
Usage Note: Setting a particular input component to zero gives the default option for that component in all cases.
The first five input components are things that every user might wish to modify; the remainder will usually be of interest only to more experienced users.
Optional every-user input:
input(1) = nnimax = maximum number of nonlinear iterations (default 200).
input(2) = ijacv = flag for determining the method of J*v evaluation:
0 => finite-difference evaluation (default)
1 => analytic evaluation
input(3) = ikrysl = flag for determining the Krylov solver:
0 => GMRES (default)
1 => BiCGSTAB
2 => TFQMR
For brief descriptions of the solvers plus references,
see the subroutines nitgm, nitstb, and nittfq.
input(4) = kdmax = maximum Krylov subspace dimension when GMRES is used
(default 20).
input(5) = irpre = flag for right preconditioning:
0 => no right preconditioning
1 => right preconditioning
Optional experienced user input:
input(6) = iksmax = maximum allowable number of iterations per call
to the Krylov solver routine (default 1000).
input(7) = iresup = residual update flag when GMRES is used; on
restarts, the residual is updated as follows:
0 => linear combination (default)
1 => direct evaluation
The first is cheap (one n-vector saxpy) but may lose
accuracy with extreme residual reduction; the second
retains accuracy better but costs one J*v product per
restart.
input(8) = ifdord = order of the finite-difference formula (sometimes)
used when input(2) = ijacv = 0. When input(2) = ijacv = 0,
this must be 0, 1, 2, or 4 on input; otherwise, it is
irrelevant. With input(2) = ijacv = 0, the precise
meaning is as follows:
If GMRES is used, then ifdord matters only if input(7) =
iresup = 1, in which case it determines the order of
the finite-difference formula used in evaluating the
initial residual at each GMRES restart (default 2); if
ifdord = 0 on input, then it is set to 2 below. NOTE: This
only affects initial residuals at restarts; first-order
differences are always used within each GMRES cycle. Using
higher-order differences at restarts only should give
the same accuracy as if higher-order differences were
used throughout; see K. Turner and H. F. Walker, "Efficient
high accuracy solutions with GMRES(m)," SIAM J. Sci.
Stat. Comput., 13 (1992), pp. 815--825.
If BiCGSTAB or TFQMR is used, then ifdord determines the
order of the finite-difference formula used at each
iteration (default 1); if ifdord = 0 on input, then it
is set to 1 below.
input(9) = ibtmax = maximum allowable number of backtracks (step
reductions) per call to nitbt (default 10).
USAGE NOTE: Backtracking can be turned off by setting
ibtmax = -1. Other negative values of ibtmax are not
valid.
input(10) = ieta = flag determining the forcing term eta as follows:
0 => abs( ||fcur|| - ||fprev+Jprev*sprev|| )/||fprev||
(default)
1 => (||fcur||/||fprev||)**2
2 => gamma*(||fcur||/||fprev||)**alpha
for user-supplied gamma in (0,1] and alpha in (1,2]
3 => fixed (constant) eta in (0,1), either 0.1 (default)
or specified by the user (see USAGE NOTE below)
Here, fcur = current f, fprev = previous f, etc. The Krylov
iterations are terminated when an iterate s satisfies
an inexact Newton condition ||F + J*s|| .le. eta*||F||.
USAGE NOTE: If input(10) = ieta = 2, then alpha and gamma
must be set in common block nitparam.h as described below.
If input(10) = ieta = 3, then the desired constant eta may
be similarly set in nitparam.h if a value other than the
default of 0.1 is desired.
The first three expressions above are from S. C. Eisenstat
and H. F. Walker, "Choosing the forcing terms in an inexact
Newton method", SIAM J. Scientific Computing, 17 (1996),
pp. 16--32. (They may be modified according to certain
safeguards in subroutine nitdrv.) The first gives convergence
that is q-superlinear and of r-order (1+sqrt(5))/2. The
second gives convergence that is r-quadratic and of q-order
p for every p in [1,2). The third gives convergence that is
of q-order alpha when gamma < 1 and, when gamma = 1, of
r-order alpha and q-order p for every p in [1,alpha). The
fourth gives q-linear convergence with asymptotic rate
constant eta in a certain norm; see R. S. Dembo, S. C.
Eisenstat, and T. Steihaug, "Inexact Newton methods",
SIAM J. Numer. Anal., 18 (1982), pp. 400-408.
Of these four choices, the 1st is usually satisfactory,
the 2nd or 3rd is sometimes preferred, and the 4th may be
useful in some situations, e.g., it may be desirable to
choose a fairly large fixed eta in (0,1), such as eta = .1,
when numerical inaccuracy prevents the Krylov solver
from obtaining much residual reduction.
On output, the components of the info argument are as follows:
info(1) = nfe (number of function evaluations) info(2) = njve (number of J*v evaluations) info(3) = nrpre (number of P(inverse)*v evaluations) info(4) = nli (number of linear iterations) info(5) = nni (number of nonlinear iterations) info(6) = nbt (number of backtracks)