24 use,
intrinsic :: iso_fortran_env
34 integer,
public,
parameter :: &
35 SK_CG = 1, & !< Conjugate Gradient Method
62 integer :: solver_type
63 integer :: krylov_size
64 integer :: preconditioning
68 real(real64) :: residual_norm
69 real(real64) :: rel_tolerance
70 real(real64) :: abs_tolerance
72 real(real64),
allocatable :: sk_work(:), sk_b(:), sk_y(:)
75 real(real64) :: fpar(16)
80 subroutine cg(n, rhs, sol, ipar, fpar, w)
81 use,
intrinsic :: iso_fortran_env
83 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
86 subroutine cgnr(n, rhs, sol, ipar, fpar, wk)
87 use,
intrinsic :: iso_fortran_env
89 real(real64) rhs(n),sol(n),fpar(16),wk(n,*)
92 subroutine bcg(n, rhs, sol, ipar, fpar, w)
93 use,
intrinsic :: iso_fortran_env
95 real(real64) fpar(16), rhs(n), sol(n), w(n,*)
98 subroutine dbcg(n, rhs, sol, ipar, fpar, w)
99 use,
intrinsic :: iso_fortran_env
101 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
104 subroutine bcgstab(n, rhs, sol, ipar, fpar, w)
105 use,
intrinsic :: iso_fortran_env
107 real(real64) rhs(n), sol(n), fpar(16), w(n,8)
110 subroutine tfqmr(n, rhs, sol, ipar, fpar, w)
111 use,
intrinsic :: iso_fortran_env
113 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
116 subroutine fom(n, rhs, sol, ipar, fpar, w)
117 use,
intrinsic :: iso_fortran_env
119 real(real64) rhs(n), sol(n), fpar(16), w(*)
122 subroutine gmres(n, rhs, sol, ipar, fpar, w)
123 use,
intrinsic :: iso_fortran_env
125 real(real64) rhs(n), sol(n), fpar(16), w(*)
128 subroutine dqgmres(n, rhs, sol, ipar, fpar, w)
129 use,
intrinsic :: iso_fortran_env
131 real(real64) rhs(n), sol(n), fpar(16), w(*)
134 subroutine fgmres(n, rhs, sol, ipar, fpar, w)
135 use,
intrinsic :: iso_fortran_env
137 real(real64) rhs(n), sol(n), fpar(16), w(*)
146 integer,
intent(in) :: n
148 logical,
intent(in) :: is_complex
150 integer :: workspace_size, m
154 sk%is_complex = is_complex
198 call parse_variable(namespace,
'SPARSKITKrylovSubspaceSize', 15, sk%krylov_size)
201 sk%preconditioning = 0
211 call parse_variable(namespace,
'SPARSKITMaxIter', 5000, sk%maxiter)
232 call parse_variable(namespace,
'SPARSKITRelTolerance', 1e-8_real64, sk%rel_tolerance)
243 call parse_variable(namespace,
'SPARSKITAbsTolerance', 1e-10_real64, sk%abs_tolerance)
252 call parse_variable(namespace,
'SPARSKITVerboseSolver', .false., sk%verbose)
266 if (mod(m, 2) /= 0) m = m + 1
268 select case (sk%solver_type)
270 message(1) =
'Info: SPARSKIT solver type: Conjugate Gradient Method'
271 workspace_size = 5*sk%size
273 message(1) =
'Info: SPARSKIT solver type: Conjugate Gradient Method (Normal Residual equation)'
274 workspace_size = 5*sk%size
276 message(1) =
'Info: SPARSKIT solver type: Bi-Conjugate Gradient Method'
277 workspace_size = 7*sk%size
279 message(1) =
'Info: SPARSKIT solver type: BCG with partial pivoting'
280 workspace_size = 11*sk%size
282 message(1) =
'Info: SPARSKIT solver type: BCG stabilized'
283 workspace_size = 8*sk%size
285 message(1) =
'Info: SPARSKIT solver type: Transpose-Free Quasi-Minimum Residual method'
286 workspace_size = 11*sk%size
288 message(1) =
'Info: SPARSKIT solver type: Full Orthogonalization Method'
289 workspace_size = (sk%size+3)*(m+2) + (m+1)*m/2
291 message(1) =
'Info: SPARSKIT solver type: Generalized Minimum Residual method'
292 workspace_size = (sk%size+3)*(m+2) + (m+1)*m/2
294 message(1) =
'Info: SPARSKIT solver type: Flexible version of Generalized Minimum Residual method'
295 workspace_size = 2*sk%size*(m+1) + (m+1)*m/2 + 3*m + 2
297 message(1) =
'Info: SPARSKIT solver type: Direct versions of Quasi-Generalized Minimum Residual method'
298 workspace_size = sk%size + (m+1) * (2*sk%size+4)
300 write(
message(1),
'(a,i4,a)')
"Input: '", sk%solver_type, &
301 "' is not a valid SPARSKIT Solver"
302 message(2) =
'( SPARSKIT Solver = cg | cgnr | bcg | dbcg | bcgstab | tfqmr | fom | gmres | fgmres | dqgmres )'
317 sk%ipar(4) = workspace_size
318 sk%ipar(5) = sk%krylov_size
321 sk%ipar(6) = sk%maxiter
324 sk%fpar(1) = sk%rel_tolerance
327 sk%fpar(2) = sk%abs_tolerance
330 safe_allocate(sk%sk_b(1:sk%size))
331 safe_allocate(sk%sk_y(1:sk%size))
332 safe_allocate(sk%sk_work(1:workspace_size))
346 safe_deallocate_a(sk%sk_b)
347 safe_deallocate_a(sk%sk_y)
348 safe_deallocate_a(sk%sk_work)
360 sko%is_complex = ski%is_complex
362 sko%solver_type = ski%solver_type
363 sko%krylov_size = ski%krylov_size
364 sko%preconditioning = ski%preconditioning
365 sko%maxiter = ski%maxiter
366 sko%used_iter = ski%used_iter
367 sko%iter_out = ski%iter_out
368 sko%residual_norm = ski%residual_norm
369 sko%rel_tolerance = ski%rel_tolerance
370 sko%abs_tolerance = ski%abs_tolerance
371 safe_allocate_source_a(sko%sk_work, ski%sk_work)
372 safe_allocate_source_a(sko%sk_b, ski%sk_b)
373 safe_allocate_source_a(sko%sk_y, ski%sk_y)
376 sko%verbose = ski%verbose
383#include "sparskit_inc.F90"
386#include "complex.F90"
387#include "sparskit_inc.F90"
real(real64), parameter, public m_zero
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine, public sparskit_solver_init(namespace, n, sk, is_complex)
integer, parameter, public sk_gmres
Generalized Minimum Residual method.
subroutine, public dsparskit_solver_run(namespace, sk, op, opt, sol, rhs)
subroutine, public zsparskit_solver_run(namespace, sk, op, opt, sol, rhs)
subroutine, public sparskit_solver_end(sk)
subroutine, public sparskit_solver_copy(sko, ski)
integer, parameter, public sk_dbcg
BCG with partial pivoting.
integer, parameter, public sk_bcg
Bi-Conjugate Gradient Method.
integer, parameter, public sk_dqgmres
Direct versions of Quasi Generalized Minimum Residual method.
integer, parameter, public sk_minval
integer, parameter, public sk_fom
Full Orthogonalization Method.
integer, parameter, public sk_bcgstab
BCG stabilized.
integer, parameter, public sk_cgnr
Conjugate Gradient Method (Normal Residual equation)
integer, parameter, public sk_maxval
integer, parameter, public sk_tfqmr
Transpose-Free Quasi-Minimum Residual method.
integer, parameter, public sk_fgmres
Flexible version of Generalized Minimum Residual method.