Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
gmres.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
34module gmres
35 use krylov, only : ksp_t, ksp_monitor_t
36 use precon, only : pc_t
37 use ax_product, only : ax_t
38 use num_types, only: rp, xp
39 use field, only : field_t
40 use coefs, only : coef_t
41 use gather_scatter, only : gs_t, gs_op_add
42 use bc_list, only : bc_list_t
43 use math, only : glsc3, rzero, rone, copy, sub2, cmult2, abscmp
44 use neko_config, only : neko_blk_size
45 use comm
46 implicit none
47 private
48
50 type, public, extends(ksp_t) :: gmres_t
51 integer :: lgmres = 30
52 real(kind=rp), allocatable :: w(:)
53 real(kind=rp), allocatable :: r(:)
54 real(kind=rp), allocatable :: z(:,:)
55 real(kind=rp), allocatable :: v(:,:)
57 real(kind=xp), allocatable :: h(:,:)
58 real(kind=xp), allocatable :: s(:)
59 real(kind=xp), allocatable :: gam(:)
60 real(kind=xp), allocatable :: c(:)
61 contains
62 procedure, pass(this) :: init => gmres_init
63 procedure, pass(this) :: free => gmres_free
64 procedure, pass(this) :: solve => gmres_solve
65 procedure, pass(this) :: solve_coupled => gmres_solve_coupled
66 end type gmres_t
67
68contains
69
71 subroutine gmres_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
72 class(gmres_t), target, intent(inout) :: this
73 integer, intent(in) :: n
74 integer, intent(in) :: max_iter
75 class(pc_t), optional, intent(in), target :: M
76 real(kind=rp), optional, intent(in) :: rel_tol
77 real(kind=rp), optional, intent(in) :: abs_tol
78 logical, optional, intent(in) :: monitor
79
80 call this%free()
81
82 if (present(m)) then
83 this%M => m
84 end if
85
86 allocate(this%w(n))
87 allocate(this%r(n))
88
89 allocate(this%c(this%lgmres))
90 allocate(this%s(this%lgmres))
91 allocate(this%gam(this%lgmres + 1))
92
93 allocate(this%z(n, this%lgmres))
94 allocate(this%v(n, this%lgmres))
95
96 allocate(this%h(this%lgmres, this%lgmres))
97
98 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
99 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
100 else if (present(rel_tol) .and. present(abs_tol)) then
101 call this%ksp_init(max_iter, rel_tol, abs_tol)
102 else if (present(monitor) .and. present(abs_tol)) then
103 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
104 else if (present(rel_tol) .and. present(monitor)) then
105 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
106 else if (present(rel_tol)) then
107 call this%ksp_init(max_iter, rel_tol = rel_tol)
108 else if (present(abs_tol)) then
109 call this%ksp_init(max_iter, abs_tol = abs_tol)
110 else if (present(monitor)) then
111 call this%ksp_init(max_iter, monitor = monitor)
112 else
113 call this%ksp_init(max_iter)
114 end if
115
116 end subroutine gmres_init
117
119 subroutine gmres_free(this)
120 class(gmres_t), intent(inout) :: this
121
122 call this%ksp_free()
123
124 if (allocated(this%w)) then
125 deallocate(this%w)
126 end if
127
128 if (allocated(this%c)) then
129 deallocate(this%c)
130 end if
131
132 if (allocated(this%r)) then
133 deallocate(this%r)
134 end if
135
136 if (allocated(this%z)) then
137 deallocate(this%z)
138 end if
139
140 if (allocated(this%h)) then
141 deallocate(this%h)
142 end if
143
144 if (allocated(this%v)) then
145 deallocate(this%v)
146 end if
147
148 if (allocated(this%s)) then
149 deallocate(this%s)
150 end if
151
152
153 if (allocated(this%gam)) then
154 deallocate(this%gam)
155 end if
156
157 nullify(this%M)
158
159 end subroutine gmres_free
160
162 function gmres_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) &
163 result(ksp_results)
164 class(gmres_t), intent(inout) :: this
165 class(ax_t), intent(in) :: ax
166 type(field_t), intent(inout) :: x
167 integer, intent(in) :: n
168 real(kind=rp), dimension(n), intent(in) :: f
169 type(coef_t), intent(inout) :: coef
170 type(bc_list_t), intent(inout) :: blst
171 type(gs_t), intent(inout) :: gs_h
172 type(ksp_monitor_t) :: ksp_results
173 integer, optional, intent(in) :: niter
174 integer :: iter, max_iter
175 integer :: i, j, k, l, ierr
176 real(kind=xp) :: w_plus(neko_blk_size), x_plus(neko_blk_size)
177 real(kind=xp) :: alpha, lr, alpha2, norm_fac
178 real(kind=rp) :: temp, rnorm
179 logical :: conv
180
181 conv = .false.
182 iter = 0
183
184 if (present(niter)) then
185 max_iter = niter
186 else
187 max_iter = this%max_iter
188 end if
189
190 associate(w => this%w, c => this%c, r => this%r, z => this%z, h => this%h, &
191 v => this%v, s => this%s, gam => this%gam)
192
193 norm_fac = 1.0_rp / sqrt(coef%volume)
194 call rzero(x%x, n)
195 gam = 0.0_xp
196 s = 1.0_xp
197 c = 1.0_xp
198 h = 0.0_xp
199 call this%monitor_start('GMRES')
200 do while (.not. conv .and. iter .lt. max_iter)
201
202 if (iter .eq. 0) then
203 call copy(r, f, n)
204 else
205 call copy(r, f, n)
206 call ax%compute(w, x%x, coef, x%msh, x%Xh)
207 call gs_h%op(w, n, gs_op_add)
208 call blst%apply(w, n)
209 call sub2(r, w, n)
210 end if
211
212 gam(1) = sqrt(glsc3(r, r, coef%mult, n))
213 if (iter .eq. 0) then
214 ksp_results%res_start = gam(1) * norm_fac
215 end if
216
217 if (abscmp(gam(1), 0.0_xp)) return
218
219 rnorm = 0.0_rp
220 temp = 1.0_rp / gam(1)
221 call cmult2(v(1,1), r, temp, n)
222 do j = 1, this%lgmres
223 iter = iter+1
224
225 call this%M%solve(z(1,j), v(1,j), n)
226
227 call ax%compute(w, z(1,j), coef, x%msh, x%Xh)
228 call gs_h%op(w, n, gs_op_add)
229 call blst%apply(w, n)
230
231 do l = 1, j
232 h(l,j) = 0.0_rp
233 end do
234
235 do i = 0, n, neko_blk_size
236 if (i + neko_blk_size .le. n) then
237 do l = 1, j
238 do k = 1, neko_blk_size
239 h(l,j) = h(l,j) + &
240 w(i+k) * v(i+k,l) * coef%mult(i+k,1,1,1)
241 end do
242 end do
243 else
244 do k = 1, n-i
245 do l = 1, j
246 h(l,j) = h(l,j) + &
247 w(i+k) * v(i+k,l) * coef%mult(i+k,1,1,1)
248 end do
249 end do
250 end if
251 end do
252
253 call mpi_allreduce(mpi_in_place, h(1,j), j, &
254 mpi_extra_precision, mpi_sum, neko_comm, ierr)
255
256 alpha2 = 0.0_rp
257 do i = 0, n, neko_blk_size
258 if (i + neko_blk_size .le. n) then
259 do k = 1, neko_blk_size
260 w_plus(k) = 0.0_rp
261 end do
262 do l = 1,j
263 do k = 1, neko_blk_size
264 w_plus(k) = w_plus(k) - h(l,j) * v(i+k,l)
265 end do
266 end do
267 do k = 1, neko_blk_size
268 w(i+k) = w(i+k) + w_plus(k)
269 alpha2 = alpha2 + w(i+k)**2 * coef%mult(i+k,1,1,1)
270 end do
271 else
272 do k = 1, n-i
273 w_plus(1) = 0.0_rp
274 do l = 1, j
275 w_plus(1) = w_plus(1) - h(l,j) * v(i+k,l)
276 end do
277 w(i+k) = w(i+k) + w_plus(1)
278 alpha2 = alpha2 + (w(i+k)**2) * coef%mult(i+k,1,1,1)
279 end do
280 end if
281 end do
282
283 call mpi_allreduce(mpi_in_place,alpha2, 1, &
284 mpi_extra_precision, mpi_sum, neko_comm, ierr)
285 alpha = sqrt(alpha2)
286 do i = 1, j-1
287 temp = h(i,j)
288 h(i,j) = c(i)*temp + s(i) * h(i+1,j)
289 h(i+1,j) = -s(i)*temp + c(i) * h(i+1,j)
290 end do
291
292 rnorm = 0.0_rp
293 if (abscmp(alpha, 0.0_xp)) then
294 conv = .true.
295 exit
296 end if
297
298 lr = sqrt(h(j,j) * h(j,j) + alpha2)
299 temp = 1.0_rp / lr
300 c(j) = h(j,j) * temp
301 s(j) = alpha * temp
302 h(j,j) = lr
303 gam(j+1) = -s(j) * gam(j)
304 gam(j) = c(j) * gam(j)
305 rnorm = abs(gam(j+1)) * norm_fac
306 call this%monitor_iter(iter, rnorm)
307 if (rnorm .lt. this%abs_tol) then
308 conv = .true.
309 exit
310 end if
311
312 if (iter + 1 .gt. max_iter) exit
313
314 if (j .lt. this%lgmres) then
315 temp = 1.0_rp / alpha
316 call cmult2(v(1,j+1), w, temp, n)
317 end if
318
319 end do
320
321 j = min(j, this%lgmres)
322 do k = j, 1, -1
323 temp = gam(k)
324 do i = j, k+1, -1
325 temp = temp - h(k,i) * c(i)
326 end do
327 c(k) = temp / h(k,k)
328 end do
329
330 do i = 0, n, neko_blk_size
331 if (i + neko_blk_size .le. n) then
332 do k = 1, neko_blk_size
333 x_plus(k) = 0.0_rp
334 end do
335 do l = 1,j
336 do k = 1, neko_blk_size
337 x_plus(k) = x_plus(k) + c(l) * z(i+k,l)
338 end do
339 end do
340 do k = 1, neko_blk_size
341 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(k)
342 end do
343 else
344 do k = 1, n-i
345 x_plus(1) = 0.0_rp
346 do l = 1, j
347 x_plus(1) = x_plus(1) + c(l) * z(i+k,l)
348 end do
349 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(1)
350 end do
351 end if
352 end do
353 end do
354
355 end associate
356 call this%monitor_stop()
357 ksp_results%res_final = rnorm
358 ksp_results%iter = iter
359 ksp_results%converged = this%is_converged(iter, rnorm)
360
361 end function gmres_solve
362
364 function gmres_solve_coupled(this, Ax, x, y, z, fx, fy, fz, &
365 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
366 class(gmres_t), intent(inout) :: this
367 class(ax_t), intent(in) :: ax
368 type(field_t), intent(inout) :: x
369 type(field_t), intent(inout) :: y
370 type(field_t), intent(inout) :: z
371 integer, intent(in) :: n
372 real(kind=rp), dimension(n), intent(in) :: fx
373 real(kind=rp), dimension(n), intent(in) :: fy
374 real(kind=rp), dimension(n), intent(in) :: fz
375 type(coef_t), intent(inout) :: coef
376 type(bc_list_t), intent(inout) :: blstx
377 type(bc_list_t), intent(inout) :: blsty
378 type(bc_list_t), intent(inout) :: blstz
379 type(gs_t), intent(inout) :: gs_h
380 type(ksp_monitor_t), dimension(3) :: ksp_results
381 integer, optional, intent(in) :: niter
382
383 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
384 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
385 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
386
387 end function gmres_solve_coupled
388
389end module gmres
390
391
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:38
type(mpi_datatype) mpi_extra_precision
Definition comm.F90:47
Defines a field.
Definition field.f90:34
Gather-scatter.
Defines various GMRES methods.
Definition gmres.f90:34
type(ksp_monitor_t) function gmres_solve(this, ax, x, f, n, coef, blst, gs_h, niter)
Standard GMRES solve.
Definition gmres.f90:164
type(ksp_monitor_t) function, dimension(3) gmres_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Standard GMRES coupled solve.
Definition gmres.f90:366
subroutine gmres_free(this)
Deallocate a standard GMRES solver.
Definition gmres.f90:120
subroutine gmres_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a standard GMRES solver.
Definition gmres.f90:72
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Definition math.f90:60
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:700
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:894
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:227
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:238
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:628
Build configurations.
integer, parameter neko_blk_size
integer, parameter, public xp
Definition num_types.f90:14
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Krylov preconditioner.
Definition precon.f90:34
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:47
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Standard preconditioned generalized minimal residual method.
Definition gmres.f90:50
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:68
Defines a canonical Krylov preconditioner.
Definition precon.f90:40