nittfq2 Interface

interface
public subroutine nittfq2(n, xcur, fcur, step, eta, f, jacv, rpar, ipar, ijacv, irpre, iksmax, ifdord, nfe, njve, nrpre, nli, r, rcgs, rtil, d, p, q, u, v, y, rwork1, rwork2, rsnrm, dinpr, dnorm, itrmks)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: n

Dimension of the problem

real(kind=r8), intent(in), dimension(n):: xcur

Array of length n containing the current value

real(kind=r8), intent(in), dimension(n):: fcur

Array of length n containing current approximate solution

real(kind=r8), intent(inout), dimension(n):: step

Vector of of length n containing trial step

real(kind=r8), intent(in) :: eta

Relative residual reduction factor

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(inout), dimension(*):: rpar

Parameter/work array passed to the f and jacv routines

integer, intent(inout), dimension(*):: ipar

Parameter/work array passed to the f and jacv routines

integer, intent(in) :: ijacv

Flag for determining method of evaluation. 0 (default) indicates finite-difference evaluation, while 1 indicates analytic evaluation.

integer, intent(in) :: irpre

Flag for right preconditioning. 0 indicates no preconditioning, while 1 inidcates right preconditioning.

integer, intent(in) :: iksmax

Maximum allowable number of TFQMR iterations

integer, intent(in) :: ifdord

Order of the finite-difference formula used in TFQMR when J*v products are evaluated using finite-differences. When ijacv = 0 on input to nitsol, ifdord is set to 1, 2, or 4 in nitsol; otherwise, it is irrelevant. When ijacv = 0 on input to this subroutine, ifdord determines the order of the finite-difference formula used at each TFQMR iteration (default 1). In this case, ijacv is set to -1 below to signal to nitjv that the order of the finite-difference formula is to be determined by ifdord. The original value ijacv = 0 is restored on return.

integer, intent(inout) :: nfe

Number of function evaluations

integer, intent(inout) :: njve

Number of evaluations

integer, intent(inout) :: nrpre

Number of evaluations

integer, intent(inout) :: nli

Number of linear iterations

real(kind=r8), intent(out), dimension(n):: r

Residual vector (for the QMR process)

real(kind=r8), intent(out), dimension(n):: rcgs

Residual vector (of the underlying CGS process)

real(kind=r8), intent(out), dimension(n):: rtil

'Shadow' residual vector used in bi-orthogonalization

real(kind=r8), intent(out), dimension(n):: d

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: p

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: q

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: u

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: v

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: y

Vector used in TFQMR

real(kind=r8), intent(out), dimension(n):: rwork1

Work vector, passed on to nitjv

real(kind=r8), intent(out), dimension(n):: rwork2

Work vector, passed on to nitjv

real(kind=r8), intent(out) :: rsnrm

TFQMR residual norm on return

procedure(dinpr_intr) :: dinpr

Inner-product routine, either user-supplied or BLAS ddot.

procedure(dnorm_intr) :: dnorm

Norm routine, either user supplied or BLAS dnrm2.

integer, intent(out) :: itrmks

Termination flag. Values have the following meanings:

0
normal termination: acceptable step found
1
failure in nitjv
2
failure in nitjv
3
Acceptable step not found in iksmax TFQMR iterations
4
TFQMR breakdown
5
Floating point error (the underlying CGS iteration has probably blown up)

Description

An interface to my modified versino of the nitsol implementation of the transpose-free quasi-minimal residual method (TFQMR) for iteratively solving linear systems. It has been modified so that the user provides a non-zero initial guess of the solution.