72 subroutine gmres_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
73 class(
gmres_t),
target,
intent(inout) :: this
74 integer,
intent(in) :: n
75 integer,
intent(in) :: max_iter
76 class(
pc_t),
optional,
intent(in),
target :: M
77 real(kind=
rp),
optional,
intent(in) :: rel_tol
78 real(kind=
rp),
optional,
intent(in) :: abs_tol
79 logical,
optional,
intent(in) :: monitor
90 allocate(this%c(this%lgmres))
91 allocate(this%s(this%lgmres))
92 allocate(this%gam(this%lgmres + 1))
94 allocate(this%z(n, this%lgmres))
95 allocate(this%v(n, this%lgmres))
97 allocate(this%h(this%lgmres, this%lgmres))
99 if (
present(rel_tol) .and.
present(abs_tol) .and.
present(monitor))
then
100 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
101 else if (
present(rel_tol) .and.
present(abs_tol))
then
102 call this%ksp_init(max_iter, rel_tol, abs_tol)
103 else if (
present(monitor) .and.
present(abs_tol))
then
104 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
105 else if (
present(rel_tol) .and.
present(monitor))
then
106 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
107 else if (
present(rel_tol))
then
108 call this%ksp_init(max_iter, rel_tol = rel_tol)
109 else if (
present(abs_tol))
then
110 call this%ksp_init(max_iter, abs_tol = abs_tol)
111 else if (
present(monitor))
then
112 call this%ksp_init(max_iter, monitor = monitor)
114 call this%ksp_init(max_iter)
163 function gmres_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) &
165 class(
gmres_t),
intent(inout) :: this
166 class(
ax_t),
intent(in) :: ax
167 type(
field_t),
intent(inout) :: x
168 integer,
intent(in) :: n
169 real(kind=
rp),
dimension(n),
intent(in) :: f
170 type(
coef_t),
intent(inout) :: coef
172 type(
gs_t),
intent(inout) :: gs_h
174 integer,
optional,
intent(in) :: niter
175 integer :: iter, max_iter
176 integer :: i, j, k, l, ierr, blk_size
178 real(kind=
xp) :: alpha, lr, alpha2, norm_fac
179 real(kind=
xp) :: hl_priv(this%lgmres)
180 real(kind=
rp) :: temp, rnorm
187 if (
present(niter))
then
190 max_iter = this%max_iter
193 associate(w => this%w, c => this%c, r => this%r, z => this%z, h => this%h, &
194 v => this%v, s => this%s, gam => this%gam)
196 norm_fac = 1.0_rp / sqrt(coef%volume)
202 call this%monitor_start(
'GMRES')
203 do while (.not. conv .and. iter .lt. max_iter)
205 if (iter .eq. 0)
then
209 call ax%compute(w, x%x, coef, x%msh, x%Xh)
210 call gs_h%op(w, n, gs_op_add)
211 call blst%apply(w, n)
215 gam(1) = sqrt(
glsc3(r, r, coef%mult, n))
216 if (iter .eq. 0)
then
217 ksp_results%res_start = gam(1) * norm_fac
220 if (
abscmp(gam(1), 0.0_xp))
exit
223 temp = 1.0_rp / gam(1)
224 call cmult2(v(1,1), r, temp, n)
225 do j = 1, this%lgmres
228 call this%M%solve(z(1,j), v(1,j), n)
230 call ax%compute(w, z(1,j), coef, x%msh, x%Xh)
231 call gs_h%op(w, n, gs_op_add)
232 call blst%apply(w, n)
246 do concurrent(k = 1:blk_size)
247 hl_priv(l) = hl_priv(l) + &
248 w(i+k) * v(i+k,l) * coef%mult(i+k,1,1,1)
256 h(l,j) = h(l,j) + hl_priv(l)
260 call mpi_allreduce(mpi_in_place, h(1,j), j, &
267 do concurrent(k = 1:blk_size)
271 do concurrent(k = 1:blk_size)
272 w_plus(k) = w_plus(k) - h(l,j) * v(i+k,l)
276 w(i+k) = w(i+k) + w_plus(k)
277 alpha2 = alpha2 + w(i+k)**2 * coef%mult(i+k,1,1,1)
281 call mpi_allreduce(mpi_in_place,alpha2, 1, &
286 h(i,j) = c(i)*temp + s(i) * h(i+1,j)
287 h(i+1,j) = -s(i)*temp + c(i) * h(i+1,j)
291 if (
abscmp(alpha, 0.0_xp))
then
296 lr = sqrt(h(j,j) * h(j,j) + alpha2)
301 gam(j+1) = -s(j) * gam(j)
302 gam(j) = c(j) * gam(j)
303 rnorm = abs(gam(j+1)) * norm_fac
304 call this%monitor_iter(iter, rnorm)
305 if (rnorm .lt. this%abs_tol)
then
310 if (iter + 1 .gt. max_iter)
exit
312 if (j .lt. this%lgmres)
then
313 temp = 1.0_rp / alpha
314 call cmult2(v(1,j+1), w, temp, n)
319 j = min(j, this%lgmres)
323 temp = temp - h(k,i) * c(i)
331 do concurrent(k = 1:blk_size)
335 do concurrent(k = 1:blk_size)
336 x_plus(k) = x_plus(k) + c(l) * z(i+k,l)
339 do concurrent(k = 1:blk_size)
340 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(k)
347 call this%monitor_stop()
348 ksp_results%res_final = rnorm
349 ksp_results%iter = iter
350 ksp_results%converged = this%is_converged(iter, rnorm)
356 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
357 class(
gmres_t),
intent(inout) :: this
358 class(ax_t),
intent(in) :: ax
359 type(field_t),
intent(inout) :: x
360 type(field_t),
intent(inout) :: y
361 type(field_t),
intent(inout) :: z
362 integer,
intent(in) :: n
363 real(kind=rp),
dimension(n),
intent(in) :: fx
364 real(kind=rp),
dimension(n),
intent(in) :: fy
365 real(kind=rp),
dimension(n),
intent(in) :: fz
366 type(coef_t),
intent(inout) :: coef
367 type(bc_list_t),
intent(inout) :: blstx
368 type(bc_list_t),
intent(inout) :: blsty
369 type(bc_list_t),
intent(inout) :: blstz
370 type(gs_t),
intent(inout) :: gs_h
371 type(ksp_monitor_t),
dimension(3) :: ksp_results
372 integer,
optional,
intent(in) :: niter
374 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
375 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
376 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
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.