Octopus
sparskit.F90
Go to the documentation of this file.
1!! Copyright (C) 2005-2006 Heiko Appel
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module sparskit_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use, intrinsic :: iso_fortran_env
27 use parser_oct_m
29
30 implicit none
31
32 private
33
34 integer, public, parameter :: &
35 SK_CG = 1, & !< Conjugate Gradient Method
36 sk_cgnr = 2, &
37 sk_bcg = 3, &
38 sk_dbcg = 4, &
39 sk_bcgstab = 5, &
40 sk_tfqmr = 6, &
41 sk_fom = 7, &
42 sk_gmres = 8, &
43 sk_fgmres = 9, &
44 sk_dqgmres = 10, &
45 sk_minval = sk_cg, &
47
48 public :: &
50
51 public :: &
57
59 private
60 logical :: is_complex
61 integer :: size
62 integer :: solver_type
63 integer :: krylov_size
64 integer :: preconditioning
65 integer :: maxiter
66 integer :: used_iter
67 integer :: iter_out
68 real(real64) :: residual_norm
69 real(real64) :: rel_tolerance
70 real(real64) :: abs_tolerance
71
72 real(real64), allocatable :: sk_work(:), sk_b(:), sk_y(:)
73
74 integer :: ipar(16)
75 real(real64) :: fpar(16)
76 logical :: verbose
77 end type sparskit_solver_t
78
79 interface
80 subroutine cg(n, rhs, sol, ipar, fpar, w)
81 use, intrinsic :: iso_fortran_env
82 integer n, ipar(16)
83 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
84 end subroutine cg
85
86 subroutine cgnr(n, rhs, sol, ipar, fpar, wk)
87 use, intrinsic :: iso_fortran_env
88 integer n, ipar(16)
89 real(real64) rhs(n),sol(n),fpar(16),wk(n,*)
90 end subroutine cgnr
91
92 subroutine bcg(n, rhs, sol, ipar, fpar, w)
93 use, intrinsic :: iso_fortran_env
94 integer n, ipar(16)
95 real(real64) fpar(16), rhs(n), sol(n), w(n,*)
96 end subroutine bcg
97
98 subroutine dbcg(n, rhs, sol, ipar, fpar, w)
99 use, intrinsic :: iso_fortran_env
100 integer n, ipar(16)
101 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
102 end subroutine dbcg
103
104 subroutine bcgstab(n, rhs, sol, ipar, fpar, w)
105 use, intrinsic :: iso_fortran_env
106 integer n, ipar(16)
107 real(real64) rhs(n), sol(n), fpar(16), w(n,8)
108 end subroutine bcgstab
109
110 subroutine tfqmr(n, rhs, sol, ipar, fpar, w)
111 use, intrinsic :: iso_fortran_env
112 integer n, ipar(16)
113 real(real64) rhs(n), sol(n), fpar(16), w(n,*)
114 end subroutine tfqmr
115
116 subroutine fom(n, rhs, sol, ipar, fpar, w)
117 use, intrinsic :: iso_fortran_env
118 integer n, ipar(16)
119 real(real64) rhs(n), sol(n), fpar(16), w(*)
120 end subroutine fom
121
122 subroutine gmres(n, rhs, sol, ipar, fpar, w)
123 use, intrinsic :: iso_fortran_env
124 integer n, ipar(16)
125 real(real64) rhs(n), sol(n), fpar(16), w(*)
126 end subroutine gmres
128 subroutine dqgmres(n, rhs, sol, ipar, fpar, w)
129 use, intrinsic :: iso_fortran_env
130 integer n, ipar(16)
131 real(real64) rhs(n), sol(n), fpar(16), w(*)
132 end subroutine dqgmres
133
134 subroutine fgmres(n, rhs, sol, ipar, fpar, w)
135 use, intrinsic :: iso_fortran_env
136 integer n, ipar(16)
137 real(real64) rhs(n), sol(n), fpar(16), w(*)
138 end subroutine fgmres
139 end interface
140
141contains
142
143 ! ---------------------------------------------------------
144 subroutine sparskit_solver_init(namespace, n, sk, is_complex)
145 type(namespace_t), intent(in) :: namespace
146 integer, intent(in) :: n
147 type(sparskit_solver_t), intent(out) :: sk
148 logical, intent(in) :: is_complex
149
150 integer :: workspace_size, m
152 push_sub(sparskit_solver_init)
154 sk%is_complex = is_complex
155 ! there might be some incompatibilities to check of real/complex and available methods?
157 !%Variable SPARSKITSolver
158 !%Type integer
159 !%Default sk_bcg
160 !%Section Math::SPARSKIT
161 !%Description
162 !% Specifies what kind of linear solver will be used.
163 !%Option sk_cg 1
164 !% Conjugate Gradient Method
165 !%Option sk_cgnr 2
166 !% Conjugate Gradient Method (Normal Residual equation)
167 !%Option sk_bcg 3
168 !% Bi-Conjugate Gradient Method
169 !%Option sk_dbcg 4
170 !% BCG with partial pivoting
171 !%Option sk_bcgstab 5
172 !% BCG stabilized
173 !%Option sk_tfqmr 6
174 !% Transpose-Free Quasi-Minimum Residual method
175 !%Option sk_fom 7
176 !% Full Orthogonalization Method
177 !%Option sk_gmres 8
178 !% Generalized Minimum Residual method
179 !%Option sk_fgmres 9
180 !% Flexible version of Generalized Minimum Residual method
181 !%Option sk_dqgmres 10
182 !% Direct versions of the Quasi-Generalized Minimum Residual method
183 !%End
184 call parse_variable(namespace, 'SPARSKITSolver', sk_bcg, sk%solver_type)
185 if (sk%solver_type < sk_minval.or.sk%solver_type > sk_maxval) then
186 call messages_input_error(namespace, 'SPARSKITSolver')
187 end if
188
189 !%Variable SPARSKITKrylovSubspaceSize
190 !%Type integer
191 !%Default 15
192 !%Section Math::SPARSKIT
193 !%Description
194 !% Some of the SPARSKIT solvers are Krylov subspace methods.
195 !% This variable determines what size the solver will use
196 !% for the subspace.
197 !%End
198 call parse_variable(namespace, 'SPARSKITKrylovSubspaceSize', 15, sk%krylov_size)
199
200 ! preconditioner not implemented
201 sk%preconditioning = 0
202
203 !%Variable SPARSKITMaxIter
204 !%Type integer
205 !%Default 50000
206 !%Section Math::SPARSKIT
207 !%Description
208 !% This variable controls the maximum number of iteration steps that
209 !% will be performed by the (iterative) linear solver.
210 !%End
211 call parse_variable(namespace, 'SPARSKITMaxIter', 5000, sk%maxiter)
212
213 !%Variable SPARSKITIterOut
214 !%Type integer
215 !%Default -1
216 !%Section Math::SPARSKIT
217 !%Description
218 !% Determines how often status info of the solver is printed.
219 !% If <= 0, will never be printed.
220 !%End
221 call parse_variable(namespace, 'SPARSKITIterOut', -1, sk%iter_out)
222
223 !%Variable SPARSKITRelTolerance
224 !%Type float
225 !%Default 1e-8
226 !%Section Math::SPARSKIT
227 !%Description
228 !% Some SPARSKIT solvers use a relative tolerance as a stopping criterion
229 !% for the iterative solution process. This variable can be used to
230 !% specify the tolerance.
231 !%End
232 call parse_variable(namespace, 'SPARSKITRelTolerance', 1e-8_real64, sk%rel_tolerance)
233
234 !%Variable SPARSKITAbsTolerance
235 !%Type float
236 !%Default 1e-10
237 !%Section Math::SPARSKIT
238 !%Description
239 !% Some SPARSKIT solvers use an absolute tolerance as a stopping criterion
240 !% for the iterative solution process. This variable can be used to
241 !% specify the tolerance.
242 !%End
243 call parse_variable(namespace, 'SPARSKITAbsTolerance', 1e-10_real64, sk%abs_tolerance)
244
245 !%Variable SPARSKITVerboseSolver
246 !%Type logical
247 !%Default no
248 !%Section Math::SPARSKIT
249 !%Description
250 !% When set to yes, the SPARSKIT solver will write more detailed output.
251 !%End
252 call parse_variable(namespace, 'SPARSKITVerboseSolver', .false., sk%verbose)
253
254 ! size of the problem
255 if (is_complex) then
256 sk%size = 2*n
257 else
258 sk%size = n
259 end if
260
261 ! initialize workspace size
262 workspace_size = 0
263
264 ! Krylov subspace size
265 m = sk%krylov_size
266 if (mod(m, 2) /= 0) m = m + 1
267
268 select case (sk%solver_type)
269 case (sk_cg)
270 message(1) = 'Info: SPARSKIT solver type: Conjugate Gradient Method'
271 workspace_size = 5*sk%size
272 case (sk_cgnr)
273 message(1) = 'Info: SPARSKIT solver type: Conjugate Gradient Method (Normal Residual equation)'
274 workspace_size = 5*sk%size
275 case (sk_bcg)
276 message(1) = 'Info: SPARSKIT solver type: Bi-Conjugate Gradient Method'
277 workspace_size = 7*sk%size
278 case (sk_dbcg)
279 message(1) = 'Info: SPARSKIT solver type: BCG with partial pivoting'
280 workspace_size = 11*sk%size
281 case (sk_bcgstab)
282 message(1) = 'Info: SPARSKIT solver type: BCG stabilized'
283 workspace_size = 8*sk%size
284 case (sk_tfqmr)
285 message(1) = 'Info: SPARSKIT solver type: Transpose-Free Quasi-Minimum Residual method'
286 workspace_size = 11*sk%size
287 case (sk_fom)
288 message(1) = 'Info: SPARSKIT solver type: Full Orthogonalization Method'
289 workspace_size = (sk%size+3)*(m+2) + (m+1)*m/2
290 case (sk_gmres)
291 message(1) = 'Info: SPARSKIT solver type: Generalized Minimum Residual method'
292 workspace_size = (sk%size+3)*(m+2) + (m+1)*m/2
293 case (sk_fgmres)
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
296 case (sk_dqgmres)
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)
299 case default
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 )'
303 call messages_fatal(2, namespace=namespace)
304 end select
305 call messages_info(1, namespace=namespace)
306
307 ! Now we initialize the arrays for the reverse communication protocol
308 sk%ipar = 0
309 sk%fpar = 0
310
311 ! A call to the solver with ipar(1) == 0 will initialize the iterative solver.
312 sk%ipar(1) = 0
313
314 ! Stopping criteria; use convergence test scheme 2
315 sk%ipar(3) = 2
316
317 sk%ipar(4) = workspace_size
318 sk%ipar(5) = sk%krylov_size
319
320 ! Maximum number of matrix-vector multiplies
321 sk%ipar(6) = sk%maxiter
322
323 ! Relative tolerance
324 sk%fpar(1) = sk%rel_tolerance
325
326 ! Absolute tolerance
327 sk%fpar(2) = sk%abs_tolerance
328
329 ! allocate and initialize work arrays
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))
333 sk%sk_work = m_zero
334 sk%sk_y = m_zero
335 sk%sk_b = m_zero
336
337 pop_sub(sparskit_solver_init)
338 end subroutine sparskit_solver_init
339
340 ! ---------------------------------------------------------
341 subroutine sparskit_solver_end(sk)
342 type(sparskit_solver_t), intent(inout) :: sk
343
344 push_sub(sparskit_solver_end)
345
346 safe_deallocate_a(sk%sk_b)
347 safe_deallocate_a(sk%sk_y)
348 safe_deallocate_a(sk%sk_work)
349
350 pop_sub(sparskit_solver_end)
351 end subroutine sparskit_solver_end
352
353 ! ---------------------------------------------------------
354 subroutine sparskit_solver_copy(sko, ski)
355 type(sparskit_solver_t), intent(inout) :: sko
356 type(sparskit_solver_t), intent(in) :: ski
357
358 push_sub(sparskit_solver_end)
359
360 sko%is_complex = ski%is_complex
361 sko%size = ski%size
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)
374 sko%ipar = ski%ipar
375 sko%fpar = ski%fpar
376 sko%verbose = ski%verbose
377
378 pop_sub(sparskit_solver_end)
379 end subroutine sparskit_solver_copy
380
381#include "undef.F90"
382#include "real.F90"
383#include "sparskit_inc.F90"
384
385#include "undef.F90"
386#include "complex.F90"
387#include "sparskit_inc.F90"
388
389! distdot function for dot products is defined in mesh_function_oct_m
390
391end module sparskit_oct_m
392
393!! Local Variables:
394!! mode: f90
395!! coding: utf-8
396!! End:
real(real64), parameter, public m_zero
Definition: global.F90:188
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
subroutine, public sparskit_solver_init(namespace, n, sk, is_complex)
Definition: sparskit.F90:238
integer, parameter, public sk_gmres
Generalized Minimum Residual method.
Definition: sparskit.F90:127
subroutine, public dsparskit_solver_run(namespace, sk, op, opt, sol, rhs)
Definition: sparskit.F90:543
subroutine, public zsparskit_solver_run(namespace, sk, op, opt, sol, rhs)
Definition: sparskit.F90:784
subroutine, public sparskit_solver_end(sk)
Definition: sparskit.F90:435
subroutine, public sparskit_solver_copy(sko, ski)
Definition: sparskit.F90:448
integer, parameter, public sk_dbcg
BCG with partial pivoting.
Definition: sparskit.F90:127
integer, parameter, public sk_bcg
Bi-Conjugate Gradient Method.
Definition: sparskit.F90:127
integer, parameter, public sk_dqgmres
Direct versions of Quasi Generalized Minimum Residual method.
Definition: sparskit.F90:127
integer, parameter, public sk_minval
Definition: sparskit.F90:127
integer, parameter, public sk_fom
Full Orthogonalization Method.
Definition: sparskit.F90:127
integer, parameter, public sk_bcgstab
BCG stabilized.
Definition: sparskit.F90:127
integer, parameter, public sk_cgnr
Conjugate Gradient Method (Normal Residual equation)
Definition: sparskit.F90:127
integer, parameter, public sk_maxval
Definition: sparskit.F90:127
integer, parameter, public sk_tfqmr
Transpose-Free Quasi-Minimum Residual method.
Definition: sparskit.F90:127
integer, parameter, public sk_fgmres
Flexible version of Generalized Minimum Residual method.
Definition: sparskit.F90:127