73 subroutine sx_gmres_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
74 class(
sx_gmres_t),
target,
intent(inout) :: this
75 integer,
intent(in) :: n
76 integer,
intent(in) :: max_iter
77 class(
pc_t),
optional,
intent(in),
target :: M
78 real(kind=
rp),
optional,
intent(in) :: rel_tol
79 real(kind=
rp),
optional,
intent(in) :: abs_tol
80 logical,
optional,
intent(in) :: monitor
94 allocate(this%c(this%lgmres))
95 allocate(this%s(this%lgmres))
96 allocate(this%gam(this%lgmres + 1))
98 allocate(this%z(n,this%lgmres))
99 allocate(this%v(n,this%lgmres))
101 allocate(this%h(this%lgmres,this%lgmres))
104 if (
present(rel_tol) .and.
present(abs_tol) .and.
present(monitor))
then
105 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
106 else if (
present(rel_tol) .and.
present(abs_tol))
then
107 call this%ksp_init(max_iter, rel_tol, abs_tol)
108 else if (
present(monitor) .and.
present(abs_tol))
then
109 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
110 else if (
present(rel_tol) .and.
present(monitor))
then
111 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
112 else if (
present(rel_tol))
then
113 call this%ksp_init(max_iter, rel_tol = rel_tol)
114 else if (
present(abs_tol))
then
115 call this%ksp_init(max_iter, abs_tol = abs_tol)
116 else if (
present(monitor))
then
117 call this%ksp_init(max_iter, monitor = monitor)
119 call this%ksp_init(max_iter)
179 function sx_gmres_solve(this, Ax, x, f, n, coef, blst, gs_h, niter)
result(ksp_results)
181 class(
ax_t),
intent(in) :: ax
182 type(
field_t),
intent(inout) :: x
183 integer,
intent(in) :: n
184 real(kind=
rp),
dimension(n),
intent(in) :: f
185 type(
coef_t),
intent(inout) :: coef
187 type(
gs_t),
intent(inout) :: gs_h
189 integer,
optional,
intent(in) :: niter
190 integer :: iter, max_iter, glb_n
191 integer :: i, j, k, ierr
192 real(kind=
rp),
parameter :: one = 1.0
193 real(kind=
rp) :: rnorm
194 real(kind=
rp) :: alpha, temp, l
195 real(kind=
rp) :: ratio, div0, norm_fac
201 glb_n = n / x%msh%nelv * x%msh%glb_nelv
203 if (
present(niter))
then
206 max_iter = this%max_iter
209 call rone(this%ml, n)
210 call rone(this%mu, n)
211 norm_fac = one / sqrt(coef%volume)
213 call rzero(this%gam, this%lgmres + 1)
214 call rone(this%s, this%lgmres)
215 call rone(this%c, this%lgmres)
216 call rzero(this%h, this%lgmres * this%lgmres)
218 call this%monitor_start(
'GMRES')
219 do while (.not. conv .and. iter .lt. max_iter)
223 call col3(this%r,this%ml,f,n)
226 call copy (this%r,f,n)
227 call ax%compute(this%w, x%x, coef, x%msh, x%Xh)
228 call gs_h%op(this%w, n, gs_op_add)
229 call blst%apply(this%w, n)
230 call add2s2(this%r,this%w,-one,n)
231 call col2(this%r,this%ml,n)
233 this%gam(1) = sqrt(
glsc3(this%r, this%r, coef%mult, n))
235 div0 = this%gam(1) * norm_fac
236 ksp_results%res_start = div0
239 if (
abscmp(this%gam(1), 0.0_rp))
return
242 temp = one / this%gam(1)
243 call cmult2(this%v(1,1), this%r, temp, n)
244 do j = 1, this%lgmres
246 call col3(this%w, this%mu, this%v(1,j), n)
249 call this%M%solve(this%z(1,j), this%w, n)
251 call ax%compute(this%w, this%z(1,j), coef, x%msh, x%Xh)
252 call gs_h%op(this%w, n, gs_op_add)
253 call blst%apply(this%w, n)
254 call col2(this%w, this%ml, n)
259 this%h(i,j) = this%h(i,j) + &
260 this%w(k) * this%v(k,i) * coef%mult(k,1,1,1)
265 call mpi_allreduce(this%h(1,j), this%wk1, j, &
267 call copy(this%h(1,j), this%wk1, j)
271 this%w(k) = this%w(k) - this%h(i,j) * this%v(k,i)
278 this%h(i ,j) = this%c(i)*temp + this%s(i)*this%h(i+1,j)
279 this%h(i+1,j) = -this%s(i)*temp + this%c(i)*this%h(i+1,j)
282 alpha = sqrt(
glsc3(this%w, this%w, coef%mult, n))
284 if(
abscmp(alpha, 0.0_rp))
then
288 l = sqrt(this%h(j,j) * this%h(j,j) + alpha**2)
290 this%c(j) = this%h(j,j) * temp
291 this%s(j) = alpha * temp
293 this%gam(j+1) = -this%s(j) * this%gam(j)
294 this%gam(j) = this%c(j) * this%gam(j)
296 rnorm = abs(this%gam(j+1)) * norm_fac
297 call this%monitor_iter(iter, rnorm)
299 if (rnorm .lt. this%abs_tol)
then
304 if (iter + 1 .gt. max_iter)
exit
306 if( j .lt. this%lgmres)
then
308 call cmult2(this%v(1,j+1), this%w, temp, n)
311 j = min(j, this%lgmres)
316 temp = temp - this%h(k,i) * this%c(i)
318 this%c(k) = temp / this%h(k,k)
323 x%x(k,1,1,1) = x%x(k,1,1,1) + this%c(i) * this%z(k,i)
327 call this%monitor_stop()
328 ksp_results%res_final = rnorm
329 ksp_results%iter = iter
330 ksp_results%converged = this%is_converged(iter, rnorm)
335 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
337 class(
ax_t),
intent(in) :: ax
338 type(
field_t),
intent(inout) :: x
339 type(
field_t),
intent(inout) :: y
340 type(
field_t),
intent(inout) :: z
341 integer,
intent(in) :: n
342 real(kind=
rp),
dimension(n),
intent(in) :: fx
343 real(kind=
rp),
dimension(n),
intent(in) :: fy
344 real(kind=
rp),
dimension(n),
intent(in) :: fz
345 type(
coef_t),
intent(inout) :: coef
349 type(
gs_t),
intent(inout) :: gs_h
351 integer,
optional,
intent(in) :: niter
353 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
354 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
355 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
type(ksp_monitor_t) function, dimension(3) sx_gmres_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Standard GMRES coupled solve.