The psfun_d_krylov_mod contains the generic call to a Krylov subspace method for the computation of , for large and sparse.
The various Krylov methods are contained in associated submodules, the idea is to have different methods, associated to the same orthogonalization method in the same submodule
Simple polynomial method based on the Arnoldi orthogonalization procedure, the method builds a basis for the Krylov subspace
Simple polynomial method based on the Arnoldi orthogonalization procedure, the method builds a basis for the Krylov subspace
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(psfun_d_serial), | intent(inout) | :: | fun | Function object |
||
type(psb_dspmat_type), | intent(in) | :: | a | Distribute sparse matrix |
||
type(psb_desc_type), | intent(in) | :: | desc_a | Descriptor for the sparse matrix |
||
type(psb_d_vect_type), | intent(inout) | :: | y | Output vector |
||
type(psb_d_vect_type), | intent(inout) | :: | x | Input vector |
||
real(kind=psb_dpk_), | intent(in) | :: | eps | Requested tolerance |
||
integer(kind=psb_ipk_), | intent(out) | :: | info | Output flag |
||
integer(kind=psb_ipk_), | intent(in), | optional | :: | itmax | Maximum number of iteration |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | itrace | Trace for logoutput |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | istop | Stop criterion |
|
integer(kind=psb_ipk_), | intent(out), | optional | :: | iter | Number of iteration |
|
real(kind=psb_dpk_), | intent(out), | optional | :: | err | Last estimate error |
|
real(kind=psb_dpk_), | intent(out), | optional | allocatable | :: | res(:) | Vector of the residuals |
Simple polynomial method based on the Lanczos orthogonalization procedure, the method builds a basis for the Krylov subspace
Simple polynomial method based on the Lanczos orthogonalization procedure, the method builds a basis for the Krylov subspace
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(psfun_d_serial), | intent(inout) | :: | fun | Function object |
||
type(psb_dspmat_type), | intent(in) | :: | a | Distribute sparse matrix |
||
type(psb_desc_type), | intent(in) | :: | desc_a | Descriptor for the sparse matrix |
||
type(psb_d_vect_type), | intent(inout) | :: | y | Output vector |
||
type(psb_d_vect_type), | intent(inout) | :: | x | Input vector |
||
real(kind=psb_dpk_), | intent(in) | :: | eps | Requested tolerance |
||
integer(kind=psb_ipk_), | intent(out) | :: | info | Output flag |
||
integer(kind=psb_ipk_), | intent(in), | optional | :: | itmax | Maximum number of iteration |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | itrace | Trace for logoutput |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | istop | Stop criterion |
|
integer(kind=psb_ipk_), | intent(out), | optional | :: | iter | Number of iteration |
|
real(kind=psb_dpk_), | intent(out), | optional | :: | err | Last estimate error |
|
real(kind=psb_dpk_), | intent(out), | optional | allocatable | :: | res(:) | Vector of the residuals |
This interface contains the Shift-and-Invert method based on Arnoldi orthogonalization.
Shift-and-invert method based on the Arnoldi orthogonalization procedure the method builds a basis for the shifted Krylov subspace
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(psfun_d_serial), | intent(inout) | :: | fun | Function object |
||
type(psb_dspmat_type), | intent(in) | :: | a | Distribute sparse matrix |
||
character(len=*), | intent(in) | :: | kryl | Krylov method for the solution of inner systems |
||
type(amg_dprec_type), | intent(inout) | :: | prec | Preconditioner for the inner method |
||
real(kind=psb_dpk_), | intent(in) | :: | tau | Shift parameter of the method |
||
type(psb_desc_type), | intent(in) | :: | desc_a | Descriptor for the sparse matrix |
||
integer(kind=psb_ipk_), | intent(in) | :: | initmax | Maximum number of iteration (inner method) |
||
integer(kind=psb_ipk_), | intent(in) | :: | initrace | Trace for logoutput (inner method) |
||
integer(kind=psb_ipk_), | intent(in) | :: | inistop | Stop criterion (inner method) |
||
type(psb_d_vect_type), | intent(inout) | :: | y | Output vector |
||
type(psb_d_vect_type), | intent(inout) | :: | x | Input vector |
||
real(kind=psb_dpk_), | intent(in) | :: | eps | Requested tolerance |
||
integer(kind=psb_ipk_), | intent(out) | :: | info | Output flag |
||
integer(kind=psb_ipk_), | intent(in), | optional | :: | itmax | Maximum number of iteration |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | itrace | Trace for logoutput |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | istop | Stop criterion |
|
integer(kind=psb_ipk_), | intent(out), | optional | :: | iter | Number of iteration |
|
real(kind=psb_dpk_), | intent(out), | optional | :: | err | Last estimate error |
|
real(kind=psb_dpk_), | intent(out), | optional | allocatable | :: | res(:) | Vector of the residuals |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer(kind=psb_ipk_), | public | :: | irst | = | 10 | Size for subspace restarting (inner method) |
|
integer(kind=psb_ipk_), | public | :: | istopc | = | 1 | Stopping criterion (inner method) |
|
integer(kind=psb_ipk_), | public | :: | itmax | = | 200 | Maximum number of iterations (inner method) |
|
integer(kind=psb_ipk_), | public | :: | itrace | = | 0 | Trace for the solution (inner method) |
|
character(len=20), | public | :: | kmethd | = | 'CG' | Method for the solution of the linear system |
|
character(len=20), | public | :: | kname | = | 'ARNOLDI' | Name of the Krylov method |
|
type(amg_dprec_type), | public | :: | prec | Preconditioner for the inner solution method |
|||
real(kind=psb_dpk_), | public | :: | tau | = | dzero | Shift for shift-and-invert methods |
procedure, public, pass(meth) :: apply => psfun_d_parallel_apply | |
procedure, public, pass(meth) :: plot => psfun_d_plot_info | |
procedure, public, pass(meth) :: precbuild => psfun_d_prec_build | |
procedure, public, pass(meth) :: precinit => psfun_d_prec_init | |
generic, public :: set => setstring, setinteger | |
procedure, public, pass(meth) :: setinteger => psfun_d_setinteger | |
procedure, public, pass(meth) :: setstring => psfun_d_setstring |
This is the generic function for applying every implemented Krylov method. The general iteration parameters (like the number of iteration, the stop criterion to be used, and the verbosity of the trace) can be passed directly to this routine. All the constitutive parameters of the actual method, and the information relative to the function are instead contained in the meth and fun objects. The Descriptor object `desc_a' contains the properties of the parallel environment.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(psfun_d_krylov), | intent(inout) | :: | meth | Krylov method object |
||
type(psfun_d_serial), | intent(inout) | :: | fun | Function object |
||
type(psb_dspmat_type), | intent(in) | :: | a | Distribute sparse matrix |
||
type(psb_desc_type), | intent(in) | :: | desc_a | Descriptor for the sparse matrix |
||
type(psb_d_vect_type), | intent(inout) | :: | y | Output vector |
||
type(psb_d_vect_type), | intent(inout) | :: | x | Input vector |
||
real(kind=psb_dpk_), | intent(in) | :: | eps | Requested tolerance |
||
integer(kind=psb_ipk_), | intent(out) | :: | info | Output flag |
||
integer(kind=psb_ipk_), | intent(in), | optional | :: | itmax | Maximum number of iteration |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | itrace | Trace for logoutput |
|
integer(kind=psb_ipk_), | intent(in), | optional | :: | istop | Stop criterion |
|
integer(kind=psb_ipk_), | intent(out), | optional | :: | iter | Number of iteration |
|
real(kind=psb_dpk_), | intent(out), | optional | :: | err | Last estimate error |
|
real(kind=psb_dpk_), | intent(out), | optional | allocatable | :: | res(:) | Vector of the residuals |
This function plots the convergence history of the Krylov method
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(psfun_d_krylov), | intent(inout) | :: | meth | Krylov method |
||
type(psfun_d_serial), | intent(inout) | :: | fun | Function object |
||
integer(kind=psb_ipk_), | intent(in) | :: | iter | Number of iteration |
||
real(kind=psb_dpk_), | intent(in), | dimension(:) | :: | res | Residual vector |
|
integer(kind=psb_ipk_), | intent(out) | :: | info | Result of the Gnuplot call |
This function builds the AMG4PSBLAS preconditioner for the inner solve in a Rational Krylov method
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(psfun_d_krylov), | intent(inout) | :: | meth | Krylov method |
||
type(psb_dspmat_type), | intent(inout) | :: | a | Sparse matrix |
||
type(psb_desc_type), | intent(inout) | :: | desc_a | Descriptor for the sparse matrix |
||
integer(kind=psb_ipk_) | :: | info | Result of the init call |
This function performs the init of the preconditioner for the inner solve in the rational Krylov method
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(psfun_d_krylov), | intent(inout) | :: | meth | Krylov method |
||
type(psb_ctxt_type), | intent(in) | :: | ctxt | Parallel context |
||
character(len=20), | intent(in) | :: | ptype | PSBLAS/AMG4PSBLAS preconditioner |
||
integer(kind=psb_ipk_) | :: | info | Result of the init call |