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