74 subroutine cg_cpld_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
75 class(
cg_cpld_t),
target,
intent(inout) :: this
76 integer,
intent(in) :: max_iter
77 class(
pc_t),
optional,
intent(in),
target :: M
78 integer,
intent(in) :: n
79 real(kind=
rp),
optional,
intent(in) :: rel_tol
80 real(kind=
rp),
optional,
intent(in) :: abs_tol
81 logical,
optional,
intent(in) :: monitor
103 if (
present(rel_tol) .and.
present(abs_tol) .and.
present(monitor))
then
104 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
105 else if (
present(rel_tol) .and.
present(abs_tol))
then
106 call this%ksp_init(max_iter, rel_tol, abs_tol)
107 else if (
present(monitor) .and.
present(abs_tol))
then
108 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
109 else if (
present(rel_tol) .and.
present(monitor))
then
110 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
111 else if (
present(rel_tol))
then
112 call this%ksp_init(max_iter, rel_tol = rel_tol)
113 else if (
present(abs_tol))
then
114 call this%ksp_init(max_iter, abs_tol = abs_tol)
115 else if (
present(monitor))
then
116 call this%ksp_init(max_iter, monitor = monitor)
118 call this%ksp_init(max_iter)
207 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
209 class(
ax_t),
intent(in) :: ax
210 type(
field_t),
intent(inout) :: x
211 type(
field_t),
intent(inout) :: y
212 type(
field_t),
intent(inout) :: z
213 integer,
intent(in) :: n
214 real(kind=
rp),
dimension(n),
intent(in) :: fx
215 real(kind=
rp),
dimension(n),
intent(in) :: fy
216 real(kind=
rp),
dimension(n),
intent(in) :: fz
217 type(
coef_t),
intent(inout) :: coef
221 type(
gs_t),
intent(inout) :: gs_h
223 integer,
optional,
intent(in) :: niter
224 integer :: i, iter, max_iter
225 real(kind=
rp) :: rnorm, rtr, rtr0, rtz2, rtz1
226 real(kind=
rp) :: beta, pap, alpha, alphm, norm_fac
228 if (
present(niter))
then
231 max_iter = this%max_iter
233 norm_fac = 1.0_rp / coef%volume
235 associate(p1 => this%p1, p2 => this%p2, p3 => this%p3, z1 => this%z1, &
236 z2 => this%z2, z3 => this%z3, r1 => this%r1, r2 => this%r2, &
237 r3 => this%r3, tmp => this%tmp, w1 => this%w1, w2 => this%w2, &
241 do concurrent(i = 1:n)
242 x%x(i,1,1,1) = 0.0_rp
243 y%x(i,1,1,1) = 0.0_rp
244 z%x(i,1,1,1) = 0.0_rp
254 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
257 rtr =
glsc3(tmp, coef%mult, coef%binv, n)
258 rnorm = sqrt(rtr)*norm_fac
259 ksp_results%res_start = rnorm
260 ksp_results%res_final = rnorm
262 if (
abscmp(rnorm, 0.0_rp))
then
263 ksp_results%converged = .true.
267 call this%monitor_start(
'cpldCG')
268 do iter = 1, max_iter
269 call this%M%solve(z1, this%r1, n)
270 call this%M%solve(z2, this%r2, n)
271 call this%M%solve(z3, this%r3, n)
274 do concurrent(i = 1:n)
275 this%tmp(i) = z1(i) * r1(i) &
280 rtz1 =
glsc2(tmp, coef%mult, n)
283 if (iter .eq. 1) beta = 0.0_rp
284 do concurrent(i = 1:n)
285 p1(i) = p1(i) * beta + z1(i)
286 p2(i) = p2(i) * beta + z2(i)
287 p3(i) = p3(i) * beta + z3(i)
290 call ax%compute_vector(w1, w2, w3, p1, p2, p3, coef, x%msh, x%Xh)
293 call gs_h%op(w1, n, gs_op_add)
294 call gs_h%op(w2, n, gs_op_add)
295 call gs_h%op(w3, n, gs_op_add)
298 call blstx%apply_scalar(w1, n)
299 call blsty%apply_scalar(w2, n)
300 call blstz%apply_scalar(w3, n)
302 do concurrent(i = 1:n)
303 tmp(i) = w1(i) * p1(i) &
308 pap =
glsc2(tmp, coef%mult, n)
312 do concurrent(i = 1:n)
313 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * p1(i)
314 y%x(i,1,1,1) = y%x(i,1,1,1) + alpha * p2(i)
315 z%x(i,1,1,1) = z%x(i,1,1,1) + alpha * p3(i)
316 r1(i) = r1(i) + alphm * w1(i)
317 r2(i) = r2(i) + alphm * w2(i)
318 r3(i) = r3(i) + alphm * w3(i)
319 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
322 rtr =
glsc3(tmp, coef%mult, coef%binv, n)
323 if (iter .eq. 1) rtr0 = rtr
324 rnorm = sqrt(rtr) * norm_fac
325 call this%monitor_iter(iter, rnorm)
326 if (rnorm .lt. this%abs_tol)
then
331 call this%monitor_stop()
332 ksp_results%res_final = rnorm
333 ksp_results%iter = iter
334 ksp_results%converged = this%is_converged(iter, rnorm)
type(ksp_monitor_t) function, dimension(3) cg_cpld_solve(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Coupled PCG solve.