Neko  0.8.1
A portable framework for high-order spectral element flow simulations
krylov_fctry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2022, 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 !
34  use cg, only : cg_t
35  use cg_sx, only : sx_cg_t
36  use cg_device, only : cg_device_t
37  use cacg, only : cacg_t
38  use pipecg, only : pipecg_t
39  use pipecg_sx, only : sx_pipecg_t
40  use pipecg_device, only : pipecg_device_t
42  use bicgstab, only : bicgstab_t
43  use gmres, only : gmres_t
44  use gmres_sx, only : sx_gmres_t
45  use gmres_device, only : gmres_device_t
46  use num_types, only : rp
47  use krylov, only : ksp_t, ksp_monitor_t
48  use precon, only : pc_t
49  use utils, only : neko_error
50  use neko_config
51  implicit none
52 
54 
55 contains
56 
58  subroutine krylov_solver_factory(ksp, n, solver, max_iter, abstol, M)
59  class(ksp_t), allocatable, target, intent(inout) :: ksp
60  integer, intent(in), value :: n
61  character(len=*), intent(in) :: solver
62  integer, intent(in) :: max_iter
63  real(kind=rp), optional :: abstol
64  class(pc_t), optional, intent(inout), target :: m
65 
66  if (allocated(ksp)) then
67  call krylov_solver_destroy(ksp)
68  deallocate(ksp)
69  end if
70  if (trim(solver) .eq. 'cg') then
71  if (neko_bcknd_sx .eq. 1) then
72  allocate(sx_cg_t::ksp)
73  else if (neko_bcknd_device .eq. 1) then
74  allocate(cg_device_t::ksp)
75  else
76  allocate(cg_t::ksp)
77  end if
78  else if (trim(solver) .eq. 'pipecg') then
79  if (neko_bcknd_sx .eq. 1) then
80  allocate(sx_pipecg_t::ksp)
81  else if (neko_bcknd_device .eq. 1) then
82  if (neko_bcknd_opencl .eq. 1) then
83  call neko_error('PipeCG not supported for OpenCL')
84  end if
85  allocate(pipecg_device_t::ksp)
86  else
87  allocate(pipecg_t::ksp)
88  end if
89  else if (trim(solver) .eq. 'fusedcg') then
90  if (neko_bcknd_device .eq. 1) then
91  if (neko_bcknd_opencl .eq. 1) then
92  call neko_error('FusedCG not supported for OpenCL')
93  end if
94  allocate(fusedcg_device_t::ksp)
95  else
96  call neko_error('FusedCG only supported for CUDA/HIP')
97  end if
98  else if (trim(solver) .eq. 'cacg') then
99  allocate(cacg_t::ksp)
100  else if (trim(solver) .eq. 'gmres') then
101  if (neko_bcknd_sx .eq. 1) then
102  allocate(sx_gmres_t::ksp)
103  else if (neko_bcknd_device .eq. 1) then
104  allocate(gmres_device_t::ksp)
105  else
106  allocate(gmres_t::ksp)
107  end if
108  else if (trim(solver) .eq. 'bicgstab') then
109  allocate(bicgstab_t::ksp)
110  else
111  call neko_error('Unknown Krylov solver '//trim(solver))
112  end if
113 
114  if (present(abstol) .and. present(m)) then
115  select type(kp => ksp)
116  type is(cg_t)
117  call kp%init(n, max_iter, m = m, abs_tol = abstol)
118  type is(sx_cg_t)
119  call kp%init(n, max_iter, m = m, abs_tol = abstol)
120  type is(cg_device_t)
121  call kp%init(n, max_iter, m = m, abs_tol = abstol)
122  type is(pipecg_t)
123  call kp%init(n, max_iter, m = m, abs_tol = abstol)
124  type is(sx_pipecg_t)
125  call kp%init(n, max_iter, m = m, abs_tol = abstol)
126  type is(pipecg_device_t)
127  call kp%init(n, max_iter, m = m, abs_tol = abstol)
128  type is(fusedcg_device_t)
129  call kp%init(n, max_iter, m = m, abs_tol = abstol)
130  type is(cacg_t)
131  call kp%init(n, max_iter, m = m, abs_tol = abstol)
132  type is(gmres_t)
133  call kp%init(n, max_iter, m = m, abs_tol = abstol)
134  type is(sx_gmres_t)
135  call kp%init(n, max_iter, m = m, abs_tol = abstol)
136  type is(gmres_device_t)
137  call kp%init(n, max_iter, m = m, abs_tol = abstol)
138  type is(bicgstab_t)
139  call kp%init(n, max_iter, m = m, abs_tol = abstol)
140  end select
141  else if (present(abstol)) then
142  select type(kp => ksp)
143  type is(cg_t)
144  call kp%init(n, max_iter, abs_tol = abstol)
145  type is(sx_cg_t)
146  call kp%init(n, max_iter, abs_tol = abstol)
147  type is(cg_device_t)
148  call kp%init(n, max_iter, abs_tol = abstol)
149  type is(pipecg_t)
150  call kp%init(n, max_iter, abs_tol = abstol)
151  type is(sx_pipecg_t)
152  call kp%init(n, max_iter, abs_tol = abstol)
153  type is (pipecg_device_t)
154  call kp%init(n, max_iter, abs_tol = abstol)
155  type is (fusedcg_device_t)
156  call kp%init(n, max_iter, abs_tol = abstol)
157  type is(cacg_t)
158  call kp%init(n, max_iter, abs_tol = abstol)
159  type is(gmres_t)
160  call kp%init(n, max_iter, abs_tol = abstol)
161  type is(sx_gmres_t)
162  call kp%init(n, max_iter, abs_tol = abstol)
163  type is(gmres_device_t)
164  call kp%init(n, max_iter, abs_tol = abstol)
165  type is(bicgstab_t)
166  call kp%init(n, max_iter, abs_tol = abstol)
167  end select
168  else if (present(m)) then
169  select type(kp => ksp)
170  type is(cg_t)
171  call kp%init(n, max_iter, m = m)
172  type is(sx_cg_t)
173  call kp%init(n, max_iter, m = m)
174  type is(cg_device_t)
175  call kp%init(n, max_iter, m = m)
176  type is(pipecg_t)
177  call kp%init(n, max_iter, m = m)
178  type is(sx_pipecg_t)
179  call kp%init(n, max_iter, m = m)
180  type is (pipecg_device_t)
181  call kp%init(n, max_iter, m = m)
182  type is (fusedcg_device_t)
183  call kp%init(n, max_iter, m = m)
184  type is(cacg_t)
185  call kp%init(n, max_iter, m = m)
186  type is(gmres_t)
187  call kp%init(n, max_iter, m = m)
188  type is(sx_gmres_t)
189  call kp%init(n, max_iter, m = m)
190  type is(gmres_device_t)
191  call kp%init(n, max_iter, m = m)
192  type is(bicgstab_t)
193  call kp%init(n, max_iter, m = m)
194  end select
195  else
196  select type(kp => ksp)
197  type is(cg_t)
198  call kp%init(n, max_iter)
199  type is(sx_cg_t)
200  call kp%init(n, max_iter)
201  type is(cg_device_t)
202  call kp%init(n, max_iter)
203  type is(pipecg_t)
204  call kp%init(n, max_iter)
205  type is(sx_pipecg_t)
206  call kp%init(n, max_iter)
207  type is (pipecg_device_t)
208  call kp%init(n, max_iter)
209  type is (fusedcg_device_t)
210  call kp%init(n, max_iter)
211  type is(cacg_t)
212  call kp%init(n, max_iter)
213  type is(gmres_t)
214  call kp%init(n, max_iter)
215  type is(sx_gmres_t)
216  call kp%init(n, max_iter)
217  type is(gmres_device_t)
218  call kp%init(n, max_iter)
219  type is(bicgstab_t)
220  call kp%init(n, max_iter)
221  end select
222  end if
223 
224  end subroutine krylov_solver_factory
225 
227  subroutine krylov_solver_destroy(ksp)
228  class(ksp_t), allocatable, intent(inout) :: ksp
229 
230  if (allocated(ksp)) then
231  select type(kp => ksp)
232  type is(cg_t)
233  call kp%free()
234  type is(sx_cg_t)
235  call kp%free()
236  type is(cg_device_t)
237  call kp%free()
238  type is(pipecg_t)
239  call kp%free()
240  type is(sx_pipecg_t)
241  call kp%free()
242  type is (pipecg_device_t)
243  call kp%free()
244  type is (fusedcg_device_t)
245  call kp%free()
246  type is(cacg_t)
247  call kp%free()
248  type is(gmres_t)
249  call kp%free()
250  type is(sx_gmres_t)
251  call kp%free()
252  type is(gmres_device_t)
253  call kp%free()
254  type is(bicgstab_t)
255  call kp%free()
256  end select
257  end if
258 
259  end subroutine krylov_solver_destroy
260 
261 end module krylov_fctry
262 
Defines various Bi-Conjugate Gradient Stabilized methods.
Definition: bicgstab.f90:34
Defines a communication avoiding Conjugate Gradient method.
Definition: cacg.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
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
subroutine, public krylov_solver_destroy(ksp)
Destroy an interative Krylov solver.
subroutine, public krylov_solver_factory(ksp, n, solver, max_iter, abstol, M)
Initialize an iterative Krylov solver.
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_device
Definition: neko_config.f90:44
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
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
Device based preconditioned conjugate gradient method.
Definition: cg_device.f90:50
Standard preconditioned conjugate gradient method (SX version)
Definition: cg_sx.f90:48
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
Type for storing initial and final residuals in a Krylov solver.
Definition: krylov.f90:55
Base abstract type for a canonical Krylov method, solving .
Definition: krylov.f90:65
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