Provides an explicit interface to the NITSOL package. Variables held in common blocks which can be used to control NITSOL are also provided here. At some point I may produce a proper object oriented interface for it.
Also present is an interface to the GMRES solver provided by NITSOL. Both a direct interface and a wrapper are provided. The wrapper offers a somewhat more general and F90-ish presentation of the routine.
Type | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|
integer | :: | iplvl | ||||
integer | :: | ipunit |
Type | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|
real(kind=r8) | :: | avrate | ||||
real(kind=r8) | :: | fcurnrm | ||||
integer | :: | instep | ||||
integer | :: | newstep | ||||
integer | :: | krystat |
Type | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|
real(kind=r8) | :: | choice1_exp | ||||
real(kind=r8) | :: | choice2_exp | ||||
real(kind=r8) | :: | choice2_coef | ||||
real(kind=r8) | :: | eta_cutoff | ||||
real(kind=r8) | :: | etamax | ||||
real(kind=r8) | :: | thmin | ||||
real(kind=r8) | :: | thmax | ||||
real(kind=r8) | :: | etafixed |
An explicit interface to the nitsol Newton iterative nonlinear solver.
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 interface to my modified versino of the nitsol implementation of the generalised minimal residual method (GMRES) for iteratively solving linear systems. It has been modified so that the user provides a non-zero initial guess of the solution.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(in), | dimension(n) | :: | fcur | Array of length |
|
real(kind=r8), | intent(inout), | dimension(n) | :: | step | Vector of of length |
|
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 |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array passed to the |
|
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 GMRES iterations |
||
integer, | intent(in) | :: | iresup | Residual update flag. On GMRES restarts, the residual can
be updated using a linear combination ( |
||
integer, | intent(in) | :: | ifdord | Order of the finite-difference formula (sometimes) used on GMRES restarts when 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, the precise meaning is as follows: With GMRES, ifdord matters only if iresup = 1, in which case |
||
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 |
||
integer, | intent(in) | :: | kdmax | Maximum Krylov subspace dimension; default 10. |
||
integer, | intent(in) | :: | kdmaxp1 | kdmax + 1 |
||
real(kind=r8), | intent(out), | dimension(n, kdmaxp1) | :: | vv | Matrix for storage of Krylov basis in GMRES; on return, the residual vector is contained in the first column. |
|
real(kind=r8), | intent(out), | dimension(kdmax, kdmax) | :: | rr | Matrix for storage of triangular matrix in GMRES. |
|
real(kind=r8), | intent(out), | dimension(kdmax) | :: | svbig | Vector for storage of estimate of singular vector of |
|
real(kind=r8), | intent(out), | dimension(kdmax) | :: | svsml | Vector for storage of estimate of singular vector of |
|
real(kind=r8), | intent(out), | dimension(kdmax) | :: | w | Vector containing right-hand side of triangular system and least-squares residual norm in GMRES. |
|
real(kind=r8), | intent(out), | dimension(n) | :: | rwork | Work array |
|
real(kind=r8), | intent(out) | :: | rsnrm | GMRES residual norm on return |
||
procedure(dinpr_intr) | :: | dinpr | Inner-product routine, either user-supplied or BLAS |
|||
procedure(dnorm_intr) | :: | dnorm | Norm routine, either user supplied or BLAS dnrm2. |
|||
integer, | intent(out) | :: | itrmks | Termination flag. Values have the following meanings:
|
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.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(in), | dimension(n) | :: | fcur | Array of length |
|
real(kind=r8), | intent(inout), | dimension(n) | :: | step | Vector of of length |
|
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 |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array passed to the |
|
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 |
|||
procedure(dnorm_intr) | :: | dnorm | Norm routine, either user supplied or BLAS dnrm2. |
|||
integer, | intent(out) | :: | itrmks | Termination flag. Values have the following meanings:
|
An interface to my modified versino of the nitsol implementation of the biconjugate gradient stabilized method (BiCGSTAB) for iteratively solving linear systems. It has been modified so that the user provides a non-zero initial guess of the solution.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(in), | dimension(n) | :: | fcur | Array of length |
|
real(kind=r8), | intent(inout), | dimension(n) | :: | step | Vector of of length |
|
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 |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array passed to the |
|
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 BiCGSTAB iterations |
||
integer, | intent(in) | :: | ifdord | Order of the finite-difference formula used in BiCGSTAB 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 BiCGSTAB 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 |
|
real(kind=r8), | intent(out), | dimension(n) | :: | rtil | vector used in BiCGSTAB |
|
real(kind=r8), | intent(out), | dimension(n) | :: | p | Vector used in BiCGSTAB |
|
real(kind=r8), | intent(out), | dimension(n) | :: | phat | Vector used in BiCGSTAB |
|
real(kind=r8), | intent(out), | dimension(n) | :: | v | Vector used in BiCGSTAB |
|
real(kind=r8), | intent(out), | dimension(n) | :: | t | Vector used in BiCGSTAB |
|
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 | BiCGSTAB residual norm on return |
||
procedure(dinpr_intr) | :: | dinpr | Inner-product routine, either user-supplied or BLAS |
|||
procedure(dnorm_intr) | :: | dnorm | Norm routine, either user supplied or BLAS dnrm2. |
|||
integer, | intent(out) | :: | itrmks | Termination flag. Values have the following meanings:
|
An interface to the BLAS routine for calculating Euclidean
inner product. This can be passed to nitsol for the
argument dinpr
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | |||
real(kind=r8), | intent(in), | dimension(*) | :: | x | ||
integer, | intent(in) | :: | sx | |||
real(kind=r8), | intent(in), | dimension(*) | :: | y | ||
integer, | intent(in) | :: | sy |
An interface to the BLAS routine for calculating Euclidean
norm. This can be passed to nitsol for the argument
dnorm
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | |||
real(kind=r8), | intent(in), | dimension(*) | :: | x | ||
integer, | intent(in) | :: | sx |
Interface for function which calculates vector inner products. This has the same interace as the BLAS routine ddot.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | The length of the vectors |
||
real(kind=r8), | intent(in), | dimension(*) | :: | x | The first input vector |
|
integer, | intent(in) | :: | sx | The stride in memory between successive elements of |
||
real(kind=r8), | intent(in), | dimension(*) | :: | y | The second input vector |
|
integer, | intent(in) | :: | sy | The stride in memory between successive elements of |
Inner product of x
and y
Interface for function which calculates vector norms. This has the same interface as the BLAS routine dnrm2.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | The length of the array |
||
real(kind=r8), | intent(in), | dimension(*) | :: | x | The input vector |
|
integer, | intent(in) | :: | sx | The stride in memory between consecutive elements of |
The vector norm of x
Interface for operations representing the multiplication of a vector by a matrix, such as that for a linear operator or a preconditioner.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(in), | dimension(:) | :: | v | The vector to be multiplied |
|
real(kind=r8), | intent(in), | dimension(:) | :: | xcur | Array containing the current estimate of the independent variables in the linear system. This may not be needed, but is provided just in case. |
|
real(kind=r8), | intent(in), | dimension(:) | :: | rhs | Array containing the right hand side of the linear system. This may not be needed, but is provided just in case. |
|
real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | Parameter/work array |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array |
|
logical, | intent(out) | :: | success | Indicates whether operation was completed succesfully |
Result of the operation
Interface for a subroutine which evaluates the function the zero of which is sought.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(out), | dimension(n) | :: | fcur | Array of length |
|
real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | Parameter/work array |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array |
|
integer, | intent(out) | :: | itrmf | Termination flag. 0 means normal termination, 1 means failure to produce f(xcur) |
Interface for a subroutine which optionally evaluates or , where is the Jacobian of and is a right preconditioning operator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(in), | dimension(n) | :: | fcur | Array of length |
|
integer, | intent(in) | :: | ijob | Integer flag indicating which product is desired. 0 indicates . 1 indicates . |
||
real(kind=r8), | intent(in), | dimension(n) | :: | v | An array of length |
|
real(kind=r8), | intent(out), | dimension(n) | :: | z | An array of length n containing the desired product on output. |
|
real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | Parameter/work array |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array |
|
integer, | intent(out) | :: | itrmjv | Termination flag. 0 indcates normal termination, 1 indicatesfailure to prodce , and 2 indicates failure to produce |
A dummy subroutine which does not apply a preconditioner or
calculate an analytic Jacobian. This can be passed to nitsol
for the argument jacv
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | |||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | ||
real(kind=r8), | intent(in), | dimension(n) | :: | fcur | ||
integer, | intent(in) | :: | ijob | |||
real(kind=r8), | intent(in), | dimension(n) | :: | v | ||
real(kind=r8), | intent(out), | dimension(n) | :: | z | ||
real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | ||
integer, | intent(inout), | dimension(*) | :: | ipar | ||
integer, | intent(out) | :: | itrmjv |
A dummy subroutine which does not calculate the function.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | Dimension of the problem |
||
real(kind=r8), | intent(in), | dimension(n) | :: | xcur | Array of length |
|
real(kind=r8), | intent(out), | dimension(n) | :: | fcur | Array of length |
|
real(kind=r8), | intent(inout), | dimension(*) | :: | rpar | Parameter/work array |
|
integer, | intent(inout), | dimension(*) | :: | ipar | Parameter/work array |
|
integer, | intent(out) | :: | itrmf | Termination flag. 0 means normal termination, 1 means failure to produce f(xcur) |
A wraper for the nitsol implementation of the generalised minimal residual method (GMRES) for iteratively solving linear systems. This provides a more general interface not specifically intended for use with Newton iteration. It also uses Fortran 90 features to provide a more convenient call signature.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(inout), | dimension(:) | :: | solution | On input, an initial guess of the solution to the linear system. On output, the iteratively determined solution. |
|
procedure(mat_mult) | :: | lhs | The linear operator on the left hand side of the linear system. |
|||
real(kind=r8), | intent(in), | dimension(:) | :: | rhs | The right hand side of the linear system being solved |
|
real(kind=r8), | intent(out) | :: | resid_norm | GMRES residual norm on return |
||
integer, | intent(out) | :: | flag | Termination flag. Values have the following meanings: |
||
integer, | intent(out), | optional | :: | nlhs | Number of evaluations of the left hand side of the system |
|
integer, | intent(out), | optional | :: | nrpre | Number of evaluations of the right-preconditioner |
|
integer, | intent(out), | optional | :: | nli | Number of iterations |
|
real(kind=r8), | intent(in), | optional | :: | tol | The tolerance for the solution. Default is |
|
procedure(mat_mult), | optional | :: | precond | A right-preconditioner which may be used to improve convergence of the solution. |
||
real(kind=r8), | intent(inout), | optional | dimension(*) | :: | rpar | Parameter/work array passed to the |
integer, | intent(inout), | optional | dimension(*) | :: | ipar | Parameter/work array passed to the |
integer, | intent(in), | optional | :: | resid_update | Residual update flag. On GMRES restarts, the residual can
be updated using a linear combination ( |
|
integer, | intent(in), | optional | :: | iter_max | Maximum allowable number of GMRES iterations. Default is 1000. |
|
integer, | intent(in), | optional | :: | krylov_dim | Maximum Krylov subspace dimension; default 10. |
|
procedure(dinpr_intr), | optional | :: | inner_prod | Inner-product routine, either user-supplied or BLAS |
||
procedure(dnorm_intr), | optional | :: | norm | Norm routine, either user supplied or BLAS dnrm2. |
A wraper for the nitsol implementation of the transpose-free quasi-minimal residual method (TFQMR) for iteratively solving linear systems. This provides a more general interface not specifically intended for use with Newton iteration. It also uses Fortran 90 features to provide a more convenient call signature.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(inout), | dimension(:) | :: | solution | On input, an initial guess of the solution to the linear system. On output, the iteratively determined solution. |
|
procedure(mat_mult) | :: | lhs | The linear operator on the left hand side of the linear system. |
|||
real(kind=r8), | intent(in), | dimension(:) | :: | rhs | The right hand side of the linear system being solved |
|
real(kind=r8), | intent(out) | :: | resid_norm | TFQMR residual norm on return |
||
integer, | intent(out) | :: | flag | Termination flag. Values have the following meanings: |
||
integer, | intent(out), | optional | :: | nlhs | Number of evaluations of the left hand side of the system |
|
integer, | intent(out), | optional | :: | nrpre | Number of evaluations of the right-preconditioner |
|
integer, | intent(out), | optional | :: | nli | Number of iterations |
|
real(kind=r8), | intent(in), | optional | :: | tol | The tolerance for the solution. Default is |
|
procedure(mat_mult), | optional | :: | precond | A right-preconditioner which may be used to improve convergence of the solution. |
||
real(kind=r8), | intent(inout), | optional | dimension(*) | :: | rpar | Parameter/work array passed to the |
integer, | intent(inout), | optional | dimension(*) | :: | ipar | Parameter/work array passed to the |
integer, | intent(in), | optional | :: | resid_update | Residual update flag. On GMRES restarts, the residual can
be updated using a linear combination ( |
|
integer, | intent(in), | optional | :: | iter_max | Maximum allowable number of TFQMR iterations. Default is 1000. |
|
integer, | intent(in), | optional | :: | krylov_dim | Maximum Krylov subspace dimension; default 10. |
|
procedure(dinpr_intr), | optional | :: | inner_prod | Inner-product routine, either user-supplied or BLAS |
||
procedure(dnorm_intr), | optional | :: | norm | Norm routine, either user supplied or BLAS dnrm2. |
A wraper for the nitsol implementation of the biconjugate gradient stabilized method (BiCGSTAB) for iteratively solving linear systems. This provides a more general interface not specifically intended for use with Newton iteration. It also uses Fortran 90 features to provide a more convenient call signature.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(inout), | dimension(:) | :: | solution | On input, an initial guess of the solution to the linear system. On output, the iteratively determined solution. |
|
procedure(mat_mult) | :: | lhs | The linear operator on the left hand side of the linear system. |
|||
real(kind=r8), | intent(in), | dimension(:) | :: | rhs | The right hand side of the linear system being solved |
|
real(kind=r8), | intent(out) | :: | resid_norm | BiCGSTAB residual norm on return |
||
integer, | intent(out) | :: | flag | Termination flag. Values have the following meanings: |
||
integer, | intent(out), | optional | :: | nlhs | Number of evaluations of the left hand side of the system |
|
integer, | intent(out), | optional | :: | nrpre | Number of evaluations of the right-preconditioner |
|
integer, | intent(out), | optional | :: | nli | Number of iterations |
|
real(kind=r8), | intent(in), | optional | :: | tol | The tolerance for the solution. Default is |
|
procedure(mat_mult), | optional | :: | precond | A right-preconditioner which may be used to improve convergence of the solution. |
||
real(kind=r8), | intent(inout), | optional | dimension(*) | :: | rpar | Parameter/work array passed to the |
integer, | intent(inout), | optional | dimension(*) | :: | ipar | Parameter/work array passed to the |
integer, | intent(in), | optional | :: | resid_update | Residual update flag. On GMRES restarts, the residual can be
updated using a linear combination ( |
|
integer, | intent(in), | optional | :: | iter_max | Maximum allowable number of TFQMR iterations. Default is 1000. |
|
integer, | intent(in), | optional | :: | krylov_dim | Maximum Krylov subspace dimension; default 10. |
|
procedure(dinpr_intr), | optional | :: | inner_prod | Inner-product routine, either user-supplied or BLAS |
||
procedure(dnorm_intr), | optional | :: | norm | Norm routine, either user supplied or BLAS dnrm2. |