Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
33submodule(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
44 use bicgstab, only : bicgstab_t
45 use gmres, only : gmres_t
46 use cheby, only : cheby_t
48 use gmres_sx, only : sx_gmres_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(9) = [character(len=20) :: &
58 "cg", &
59 "pipecg", &
60 "fusedcg", &
61 "cacg", &
62 "gmres", &
63 "cheby", &
64 "bicgstab", &
65 "fusedcoupledcg", &
66 "coupledcg"]
67
68contains
69
78 module subroutine krylov_solver_factory(object, n, type_name, &
79 max_iter, abstol, m, monitor)
80 class(ksp_t), allocatable, intent(inout) :: object
81 integer, intent(in), value :: n
82 character(len=*), intent(in) :: type_name
83 integer, intent(in) :: max_iter
84 real(kind=rp), optional :: abstol
85 class(pc_t), optional, intent(in), target :: M
86 logical, optional, intent(in) :: monitor
87 character(len=:), allocatable :: type_string
88
89 if (allocated(object)) then
90 call krylov_solver_destroy(object)
91 deallocate(object)
92 end if
93
94 if (trim(type_name) .eq. 'cg') then
95 if (neko_bcknd_sx .eq. 1) then
96 allocate(sx_cg_t::object)
97 else if (neko_bcknd_device .eq. 1) then
98 allocate(cg_device_t::object)
99 else
100 allocate(cg_t::object)
101 end if
102 else if (trim(type_name) .eq. 'coupledcg') then
103 allocate(cg_cpld_t::object)
104 if (neko_bcknd_device .eq. 1) then
105 call neko_error('Coupled CG only supported for CPU')
106 end if
107 else if (trim(type_name) .eq. 'pipecg') then
108 if (neko_bcknd_sx .eq. 1) then
109 allocate(sx_pipecg_t::object)
110 else if (neko_bcknd_device .eq. 1) then
111 if (neko_bcknd_opencl .eq. 1) then
112 call neko_error('PipeCG not supported for OpenCL')
113 end if
114 allocate(pipecg_device_t::object)
115 else
116 allocate(pipecg_t::object)
117 end if
118 else if (trim(type_name) .eq. 'fusedcg') then
119 if (neko_bcknd_device .eq. 1) then
120 if (neko_bcknd_opencl .eq. 1) then
121 call neko_error('FusedCG not supported for OpenCL')
122 end if
123 allocate(fusedcg_device_t::object)
124 else
125 call neko_error('FusedCG only supported for CUDA/HIP')
126 end if
127 else if (trim(type_name) .eq. 'fusedcoupledcg') then
128 if (neko_bcknd_device .eq. 1) then
129 if (neko_bcknd_opencl .eq. 1) then
130 call neko_error('Coupled FusedCG not supported for OpenCL')
131 end if
132 allocate(fusedcg_cpld_device_t::object)
133 else
134 call neko_error('Coupled FusedCG only supported for CUDA/HIP')
135 end if
136 else if (trim(type_name) .eq. 'cacg') then
137 allocate(cacg_t::object)
138 else if (trim(type_name) .eq. 'gmres') then
139 if (neko_bcknd_sx .eq. 1) then
140 allocate(sx_gmres_t::object)
141 else if (neko_bcknd_device .eq. 1) then
142 allocate(gmres_device_t::object)
143 else
144 allocate(gmres_t::object)
145 end if
146 else if (trim(type_name) .eq. 'cheby') then
147 if (neko_bcknd_device .eq. 1) then
148 allocate(cheby_device_t::object)
149 else
150 allocate(cheby_t::object)
151 end if
152 else if (trim(type_name) .eq. 'bicgstab') then
153 allocate(bicgstab_t::object)
154 else
155 type_string = concat_string_array(ksp_known_types,&
156 new_line('A') // "- ", .true.)
157 call neko_error("Unknown Krylov solver type: " &
158 // trim(type_name) // ". Known types are: " &
159 // type_string)
160 end if
161
162 ! This select type is in principle not necessary,but we have it due to
163 ! issues with compilers, when it was not there. However, at some point we
164 ! should check if we can get away with just having one obj%init statement.
165 ! Same applies to the code in the "destroy" routine below.
166 select type (obj => object)
167 type is (cg_t)
168 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
169 monitor = monitor)
170 type is (sx_cg_t)
171 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
172 monitor = monitor)
173 type is (cg_cpld_t)
174 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
175 monitor = monitor)
176 type is (cg_device_t)
177 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
178 monitor = monitor)
179 type is (pipecg_t)
180 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
181 monitor = monitor)
182 type is (sx_pipecg_t)
183 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
184 monitor = monitor)
185 type is (pipecg_device_t)
186 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
187 monitor = monitor)
188 type is (fusedcg_device_t)
189 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
190 monitor = monitor)
191 type is (fusedcg_cpld_device_t)
192 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
193 monitor = monitor)
194 type is (cacg_t)
195 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
196 monitor = monitor)
197 type is (gmres_t)
198 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
199 monitor = monitor)
200 type is (sx_gmres_t)
201 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
202 monitor = monitor)
203 type is (gmres_device_t)
204 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
205 monitor = monitor)
206 type is (bicgstab_t)
207 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
208 monitor = monitor)
209 type is (cheby_t)
210 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
211 monitor = monitor)
212 type is (cheby_device_t)
213 call obj%init(n, max_iter, m = m, abs_tol = abstol,&
214 monitor = monitor)
215 end select
216
217 end subroutine krylov_solver_factory
218
220 module subroutine krylov_solver_destroy(object)
221 class(ksp_t), allocatable, intent(inout) :: object
222 call object%free()
223 end subroutine krylov_solver_destroy
224
225end submodule krylov_fctry
226
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.
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.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_opencl
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:276
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.
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