Neko  0.8.99
A portable framework for high-order spectral element flow simulations
krylov_fctry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
33 submodule(krylov) krylov_fctry
34  use cg, only : cg_t
35  use cg_sx, only : sx_cg_t
36  use cg_cpld, only : cg_cpld_t
37  use cg_device, only : cg_device_t
38  use cacg, only : cacg_t
39  use pipecg, only : pipecg_t
40  use pipecg_sx, only : sx_pipecg_t
41  use pipecg_device, only : pipecg_device_t
44  use bicgstab, only : bicgstab_t
45  use gmres, only : gmres_t
46  use cheby, only : cheby_t
47  use cheby_device, only : cheby_device_t
48  use gmres_sx, only : sx_gmres_t
49  use gmres_device, only : gmres_device_t
50  use num_types, only : rp
51  use precon, only : pc_t
52  use utils, only : concat_string_array
54  implicit none
55 
56  ! List of all possible types created by the factory routine
57  character(len=20) :: KSP_KNOWN_TYPES(8) = [character(len=20) :: &
58  "cg", &
59  "pipecg", &
60  "fusedcg", &
61  "cacg", &
62  "gmres", &
63  "cheby", &
64  "bicgstab", &
65  "cpldcg"]
66 
67 contains
68 
77  module subroutine krylov_solver_factory(object, n, type_name, &
78  max_iter, abstol, m, monitor)
79  class(ksp_t), allocatable, target, intent(inout) :: object
80  integer, intent(in), value :: n
81  character(len=*), intent(in) :: type_name
82  integer, intent(in) :: max_iter
83  real(kind=rp), optional :: abstol
84  class(pc_t), optional, intent(inout), target :: M
85  logical, optional, intent(in) :: monitor
86  character(len=:), allocatable :: type_string
87 
88  if (allocated(object)) then
89  call krylov_solver_destroy(object)
90  deallocate(object)
91  end if
92 
93  if (trim(type_name) .eq. 'cg') then
94  if (neko_bcknd_sx .eq. 1) then
95  allocate(sx_cg_t::object)
96  else if (neko_bcknd_device .eq. 1) then
97  allocate(cg_device_t::object)
98  else
99  allocate(cg_t::object)
100  end if
101  else if (trim(type_name) .eq. 'cpldcg') then
102  allocate(cg_cpld_t::object)
103  if (neko_bcknd_device .eq. 1) then
104  call neko_error('Coupled CG only supported for CPU')
105  end if
106  else if (trim(type_name) .eq. 'pipecg') then
107  if (neko_bcknd_sx .eq. 1) then
108  allocate(sx_pipecg_t::object)
109  else if (neko_bcknd_device .eq. 1) then
110  if (neko_bcknd_opencl .eq. 1) then
111  call neko_error('PipeCG not supported for OpenCL')
112  end if
113  allocate(pipecg_device_t::object)
114  else
115  allocate(pipecg_t::object)
116  end if
117  else if (trim(type_name) .eq. 'fusedcg') then
118  if (neko_bcknd_device .eq. 1) then
119  if (neko_bcknd_opencl .eq. 1) then
120  call neko_error('FusedCG not supported for OpenCL')
121  end if
122  allocate(fusedcg_device_t::object)
123  else
124  call neko_error('FusedCG only supported for CUDA/HIP')
125  end if
126  else if (trim(type_name) .eq. 'fcpldcg') then
127  if (neko_bcknd_device .eq. 1) then
128  if (neko_bcknd_opencl .eq. 1) then
129  call neko_error('Coupled FusedCG not supported for OpenCL')
130  end if
131  allocate(fusedcg_cpld_device_t::object)
132  else
133  call neko_error('Coupled FusedCG only supported for CUDA/HIP')
134  end if
135  else if (trim(type_name) .eq. 'cacg') then
136  allocate(cacg_t::object)
137  else if (trim(type_name) .eq. 'gmres') then
138  if (neko_bcknd_sx .eq. 1) then
139  allocate(sx_gmres_t::object)
140  else if (neko_bcknd_device .eq. 1) then
141  allocate(gmres_device_t::object)
142  else
143  allocate(gmres_t::object)
144  end if
145  else if (trim(type_name) .eq. 'cheby') then
146  if (neko_bcknd_device .eq. 1) then
147  allocate(cheby_device_t::object)
148  else
149  allocate(cheby_t::object)
150  end if
151  else if (trim(type_name) .eq. 'bicgstab') then
152  allocate(bicgstab_t::object)
153  else
154  type_string = concat_string_array(ksp_known_types,&
155  new_line('A') // "- ", .true.)
156  call neko_error("Unknown Krylov solver type: " &
157  // trim(type_name) // ". Known types are: " &
158  // type_string)
159  end if
160 
161  ! This select type is in principle not necessary,but we have it due to
162  ! issues with compilers, when it was not there. However, at some point we
163  ! should check if we can get away with just having one obj%init statement.
164  ! Same applies to the code in the "destroy" routine below.
165  if (present(abstol) .and. present(m) .and. present(monitor)) then
166  select type (obj => object)
167  type is (cg_t)
168  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
169  type is (sx_cg_t)
170  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
171  type is (cg_cpld_t)
172  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
173  type is (cg_device_t)
174  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
175  type is (pipecg_t)
176  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
177  type is (sx_pipecg_t)
178  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
179  type is (pipecg_device_t)
180  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
181  type is (fusedcg_device_t)
182  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
183  type is (fusedcg_cpld_device_t)
184  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
185  type is (cacg_t)
186  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
187  type is (gmres_t)
188  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
189  type is (sx_gmres_t)
190  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
191  type is (gmres_device_t)
192  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
193  type is (bicgstab_t)
194  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
195  type is (cheby_t)
196  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
197  type is (cheby_device_t)
198  call obj%init(n, max_iter, m = m, abs_tol = abstol, monitor = monitor)
199  end select
200  else if (present(abstol) .and. present(m)) then
201  select type (obj => object)
202  type is (cg_t)
203  call obj%init(n, max_iter, m = m, abs_tol = abstol)
204  type is (sx_cg_t)
205  call obj%init(n, max_iter, m = m, abs_tol = abstol)
206  type is (cg_cpld_t)
207  call obj%init(n, max_iter, m = m, abs_tol = abstol)
208  type is (cg_device_t)
209  call obj%init(n, max_iter, m = m, abs_tol = abstol)
210  type is (pipecg_t)
211  call obj%init(n, max_iter, m = m, abs_tol = abstol)
212  type is (sx_pipecg_t)
213  call obj%init(n, max_iter, m = m, abs_tol = abstol)
214  type is (pipecg_device_t)
215  call obj%init(n, max_iter, m = m, abs_tol = abstol)
216  type is (fusedcg_device_t)
217  call obj%init(n, max_iter, m = m, abs_tol = abstol)
218  type is (fusedcg_cpld_device_t)
219  call obj%init(n, max_iter, m = m, abs_tol = abstol)
220  type is (cacg_t)
221  call obj%init(n, max_iter, m = m, abs_tol = abstol)
222  type is (gmres_t)
223  call obj%init(n, max_iter, m = m, abs_tol = abstol)
224  type is (sx_gmres_t)
225  call obj%init(n, max_iter, m = m, abs_tol = abstol)
226  type is (gmres_device_t)
227  call obj%init(n, max_iter, m = m, abs_tol = abstol)
228  type is (bicgstab_t)
229  call obj%init(n, max_iter, m = m, abs_tol = abstol)
230  type is (cheby_t)
231  call obj%init(n, max_iter, m = m, abs_tol = abstol)
232  type is (cheby_device_t)
233  call obj%init(n, max_iter, m = m, abs_tol = abstol)
234  end select
235  else if (present(monitor) .and. present(m)) then
236  select type (obj => object)
237  type is (cg_t)
238  call obj%init(n, max_iter, m = m, monitor = monitor)
239  type is (sx_cg_t)
240  call obj%init(n, max_iter, m = m, monitor = monitor)
241  type is (cg_cpld_t)
242  call obj%init(n, max_iter, m = m, monitor = monitor)
243  type is (cg_device_t)
244  call obj%init(n, max_iter, m = m, monitor = monitor)
245  type is (pipecg_t)
246  call obj%init(n, max_iter, m = m, monitor = monitor)
247  type is (sx_pipecg_t)
248  call obj%init(n, max_iter, m = m, monitor = monitor)
249  type is (pipecg_device_t)
250  call obj%init(n, max_iter, m = m, monitor = monitor)
251  type is (fusedcg_device_t)
252  call obj%init(n, max_iter, m = m, monitor = monitor)
253  type is (fusedcg_cpld_device_t)
254  call obj%init(n, max_iter, m = m, monitor = monitor)
255  type is (cacg_t)
256  call obj%init(n, max_iter, m = m, monitor = monitor)
257  type is (gmres_t)
258  call obj%init(n, max_iter, m = m, monitor = monitor)
259  type is (sx_gmres_t)
260  call obj%init(n, max_iter, m = m, monitor = monitor)
261  type is (gmres_device_t)
262  call obj%init(n, max_iter, m = m, monitor = monitor)
263  type is (bicgstab_t)
264  call obj%init(n, max_iter, m = m, monitor = monitor)
265  type is (cheby_t)
266  call obj%init(n, max_iter, m = m, monitor = monitor)
267  type is (cheby_device_t)
268  call obj%init(n, max_iter, m = m, monitor = monitor)
269  end select
270  else if (present(abstol) .and. present(monitor)) then
271  select type (obj => object)
272  type is (cg_t)
273  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
274  type is (sx_cg_t)
275  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
276  type is (cg_cpld_t)
277  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
278  type is (cg_device_t)
279  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
280  type is (pipecg_t)
281  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
282  type is (sx_pipecg_t)
283  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
284  type is (pipecg_device_t)
285  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
286  type is (fusedcg_device_t)
287  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
288  type is (fusedcg_cpld_device_t)
289  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
290  type is (cacg_t)
291  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
292  type is (gmres_t)
293  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
294  type is (sx_gmres_t)
295  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
296  type is (gmres_device_t)
297  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
298  type is (bicgstab_t)
299  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
300  type is (cheby_t)
301  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
302  type is (cheby_device_t)
303  call obj%init(n, max_iter, monitor = monitor, abs_tol = abstol)
304  end select
305  else if (present(abstol)) then
306  select type (obj => object)
307  type is (cg_t)
308  call obj%init(n, max_iter, abs_tol = abstol)
309  type is (sx_cg_t)
310  call obj%init(n, max_iter, abs_tol = abstol)
311  type is (cg_cpld_t)
312  call obj%init(n, max_iter, abs_tol = abstol)
313  type is (cg_device_t)
314  call obj%init(n, max_iter, abs_tol = abstol)
315  type is (pipecg_t)
316  call obj%init(n, max_iter, abs_tol = abstol)
317  type is (sx_pipecg_t)
318  call obj%init(n, max_iter, abs_tol = abstol)
319  type is (pipecg_device_t)
320  call obj%init(n, max_iter, abs_tol = abstol)
321  type is (fusedcg_device_t)
322  call obj%init(n, max_iter, abs_tol = abstol)
323  type is (fusedcg_cpld_device_t)
324  call obj%init(n, max_iter, abs_tol = abstol)
325  type is (cacg_t)
326  call obj%init(n, max_iter, abs_tol = abstol)
327  type is (gmres_t)
328  call obj%init(n, max_iter, abs_tol = abstol)
329  type is (sx_gmres_t)
330  call obj%init(n, max_iter, abs_tol = abstol)
331  type is (gmres_device_t)
332  call obj%init(n, max_iter, abs_tol = abstol)
333  type is (bicgstab_t)
334  call obj%init(n, max_iter, abs_tol = abstol)
335  type is (cheby_t)
336  call obj%init(n, max_iter, abs_tol = abstol)
337  type is (cheby_device_t)
338  call obj%init(n, max_iter, abs_tol = abstol)
339  end select
340  else if (present(monitor)) then
341  select type (obj => object)
342  type is (cg_t)
343  call obj%init(n, max_iter, monitor = monitor)
344  type is (sx_cg_t)
345  call obj%init(n, max_iter, monitor = monitor)
346  type is (cg_cpld_t)
347  call obj%init(n, max_iter, monitor = monitor)
348  type is (cg_device_t)
349  call obj%init(n, max_iter, monitor = monitor)
350  type is (pipecg_t)
351  call obj%init(n, max_iter, monitor = monitor)
352  type is (sx_pipecg_t)
353  call obj%init(n, max_iter, monitor = monitor)
354  type is (pipecg_device_t)
355  call obj%init(n, max_iter, monitor = monitor)
356  type is (fusedcg_device_t)
357  call obj%init(n, max_iter, monitor = monitor)
358  type is (fusedcg_cpld_device_t)
359  call obj%init(n, max_iter, monitor = monitor)
360  type is (cacg_t)
361  call obj%init(n, max_iter, monitor = monitor)
362  type is (gmres_t)
363  call obj%init(n, max_iter, monitor = monitor)
364  type is (sx_gmres_t)
365  call obj%init(n, max_iter, monitor = monitor)
366  type is (gmres_device_t)
367  call obj%init(n, max_iter, monitor = monitor)
368  type is (bicgstab_t)
369  call obj%init(n, max_iter, monitor = monitor)
370  type is (cheby_t)
371  call obj%init(n, max_iter, monitor = monitor)
372  type is (cheby_device_t)
373  call obj%init(n, max_iter, monitor = monitor)
374  end select
375  else if (present(m)) then
376  select type (obj => object)
377  type is (cg_t)
378  call obj%init(n, max_iter, m = m)
379  type is (sx_cg_t)
380  call obj%init(n, max_iter, m = m)
381  type is (cg_cpld_t)
382  call obj%init(n, max_iter, m = m)
383  type is (cg_device_t)
384  call obj%init(n, max_iter, m = m)
385  type is (pipecg_t)
386  call obj%init(n, max_iter, m = m)
387  type is (sx_pipecg_t)
388  call obj%init(n, max_iter, m = m)
389  type is (pipecg_device_t)
390  call obj%init(n, max_iter, m = m)
391  type is (fusedcg_device_t)
392  call obj%init(n, max_iter, m = m)
393  type is (fusedcg_cpld_device_t)
394  call obj%init(n, max_iter, m = m)
395  type is (cacg_t)
396  call obj%init(n, max_iter, m = m)
397  type is (gmres_t)
398  call obj%init(n, max_iter, m = m)
399  type is (sx_gmres_t)
400  call obj%init(n, max_iter, m = m)
401  type is (gmres_device_t)
402  call obj%init(n, max_iter, m = m)
403  type is (bicgstab_t)
404  call obj%init(n, max_iter, m = m)
405  type is (cheby_t)
406  call obj%init(n, max_iter, m = m)
407  type is (cheby_device_t)
408  call obj%init(n, max_iter, m = m)
409  end select
410  else
411  select type (obj => object)
412  type is (cg_t)
413  call obj%init(n, max_iter)
414  type is (sx_cg_t)
415  call obj%init(n, max_iter)
416  type is (cg_cpld_t)
417  call obj%init(n, max_iter)
418  type is (cg_device_t)
419  call obj%init(n, max_iter)
420  type is (pipecg_t)
421  call obj%init(n, max_iter)
422  type is (sx_pipecg_t)
423  call obj%init(n, max_iter)
424  type is (pipecg_device_t)
425  call obj%init(n, max_iter)
426  type is (fusedcg_device_t)
427  call obj%init(n, max_iter)
428  type is (fusedcg_cpld_device_t)
429  call obj%init(n, max_iter)
430  type is (cacg_t)
431  call obj%init(n, max_iter)
432  type is (gmres_t)
433  call obj%init(n, max_iter)
434  type is (sx_gmres_t)
435  call obj%init(n, max_iter)
436  type is (gmres_device_t)
437  call obj%init(n, max_iter)
438  type is (bicgstab_t)
439  call obj%init(n, max_iter)
440  type is (cheby_t)
441  call obj%init(n, max_iter)
442  type is (cheby_device_t)
443  call obj%init(n, max_iter)
444  end select
445  end if
446 
447  end subroutine krylov_solver_factory
448 
450  module subroutine krylov_solver_destroy(object)
451  class(ksp_t), allocatable, intent(inout) :: object
452 
453  if (allocated(object)) then
454  select type (obj => object)
455  type is (cg_t)
456  call obj%free()
457  type is (sx_cg_t)
458  call obj%free()
459  type is (cg_cpld_t)
460  call obj%free()
461  type is (cg_device_t)
462  call obj%free()
463  type is (pipecg_t)
464  call obj%free()
465  type is (sx_pipecg_t)
466  call obj%free()
467  type is (pipecg_device_t)
468  call obj%free()
469  type is (fusedcg_device_t)
470  call obj%free()
471  type is (fusedcg_cpld_device_t)
472  call obj%free()
473  type is (cacg_t)
474  call obj%free()
475  type is (gmres_t)
476  call obj%free()
477  type is (sx_gmres_t)
478  call obj%free()
479  type is (gmres_device_t)
480  call obj%free()
481  type is (bicgstab_t)
482  call obj%free()
483  type is (cheby_t)
484  call obj%free()
485  end select
486  end if
487 
488  end subroutine krylov_solver_destroy
489 
490 end submodule krylov_fctry
491 
Defines various Bi-Conjugate Gradient Stabilized methods.
Definition: bicgstab.f90:34
Defines a communication avoiding Conjugate Gradient method.
Definition: cacg.f90:34
Defines a coupled Conjugate Gradient methods.
Definition: cg_coupled.f90:34
Defines various Conjugate Gradient methods for accelerators.
Definition: cg_device.f90:34
Defines various Conjugate Gradient methods.
Definition: cg_sx.f90:34
Defines various Conjugate Gradient methods.
Definition: cg.f90:34
Chebyshev preconditioner.
Chebyshev preconditioner.
Definition: cheby.f90:34
Defines a fused Conjugate Gradient method for accelerators.
Defines a fused Conjugate Gradient method for accelerators.
Defines various GMRES methods.
Defines various GMRES methods.
Definition: gmres_sx.f90:34
Defines various GMRES methods.
Definition: gmres.f90:34
Implements the base abstract type for Krylov solvers plus helper types.
Definition: krylov.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_sx
Definition: neko_config.f90:39
integer, parameter neko_bcknd_opencl
Definition: neko_config.f90:43
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a pipelined Conjugate Gradient methods.
Defines a pipelined Conjugate Gradient methods SX-Aurora backend.
Definition: pipecg_sx.f90:34
Defines a pipelined Conjugate Gradient methods.
Definition: pipecg.f90:34
Krylov preconditioner.
Definition: precon.f90:34
Utilities.
Definition: utils.f90:35
character(:) function, allocatable, public concat_string_array(array, sep, prepend)
Concatenate an array of strings into one string with array items separated by spaces.
Definition: utils.f90:208
Standard preconditioned Bi-Conjugate Gradient Stabilized method.
Definition: bicgstab.f90:50
S-step communication avoiding preconditioned conjugate gradient method.
Definition: cacg.f90:51
Standard preconditioned conjugate gradient method.
Definition: cg.f90:51
Coupled preconditioned conjugate gradient method.
Definition: cg_coupled.f90:49
Device based preconditioned conjugate gradient method.
Definition: cg_device.f90:50
Standard preconditioned conjugate gradient method (SX version)
Definition: cg_sx.f90:48
Defines a Chebyshev preconditioner.
Definition: cheby.f90:52
Defines a Chebyshev preconditioner.
Fused preconditioned conjugate gradient method.
Fused preconditioned conjugate gradient method.
Standard preconditioned generalized minimal residual method.
Definition: gmres.f90:49
Standard preconditioned generalized minimal residual method.
Standard preconditioned generalized minimal residual method (SX version)
Definition: gmres_sx.f90:49
Pipelined preconditioned conjugate gradient method.
Definition: pipecg.f90:51
Pipelined preconditioned conjugate gradient method.
Pipelined preconditioned conjugate gradient method for SX-Aurora.
Definition: pipecg_sx.f90:49
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40