73 subroutine cg_cpld_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
75 integer,
intent(in) :: max_iter
76 class(
pc_t),
optional,
intent(in),
target :: M
77 integer,
intent(in) :: n
78 real(kind=
rp),
optional,
intent(in) :: rel_tol
79 real(kind=
rp),
optional,
intent(in) :: abs_tol
80 logical,
optional,
intent(in) :: monitor
102 if (
present(rel_tol) .and.
present(abs_tol) .and.
present(monitor))
then
103 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
104 else if (
present(rel_tol) .and.
present(abs_tol))
then
105 call this%ksp_init(max_iter, rel_tol, abs_tol)
106 else if (
present(monitor) .and.
present(abs_tol))
then
107 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
108 else if (
present(rel_tol) .and.
present(monitor))
then
109 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
110 else if (
present(rel_tol))
then
111 call this%ksp_init(max_iter, rel_tol = rel_tol)
112 else if (
present(abs_tol))
then
113 call this%ksp_init(max_iter, abs_tol = abs_tol)
114 else if (
present(monitor))
then
115 call this%ksp_init(max_iter, monitor = monitor)
117 call this%ksp_init(max_iter)
206 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
208 class(
ax_t),
intent(in) :: ax
209 type(
field_t),
intent(inout) :: x
210 type(
field_t),
intent(inout) :: y
211 type(
field_t),
intent(inout) :: z
212 integer,
intent(in) :: n
213 real(kind=
rp),
dimension(n),
intent(in) :: fx
214 real(kind=
rp),
dimension(n),
intent(in) :: fy
215 real(kind=
rp),
dimension(n),
intent(in) :: fz
216 type(
coef_t),
intent(inout) :: coef
220 type(
gs_t),
intent(inout) :: gs_h
222 integer,
optional,
intent(in) :: niter
223 real(kind=
rp),
parameter :: one = 1.0
224 real(kind=
rp),
parameter :: zero = 0.0
225 integer :: i, iter, max_iter
226 real(kind=
rp) :: rnorm, rtr, rtr0, rtz2, rtz1
227 real(kind=
rp) :: beta, pap, alpha, alphm, norm_fac
229 if (
present(niter))
then
232 max_iter = this%max_iter
234 norm_fac = one / coef%volume
236 associate(p1 => this%p1, p2 => this%p2, p3 => this%p3, z1 => this%z1, &
237 z2 => this%z2, z3 => this%z3, r1 => this%r1, r2 => this%r2, &
238 r3 => this%r3, tmp => this%tmp, w1 => this%w1, w2 => this%w2, &
242 do concurrent(i = 1:n)
243 x%x(i,1,1,1) = 0.0_rp
244 y%x(i,1,1,1) = 0.0_rp
245 z%x(i,1,1,1) = 0.0_rp
255 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
258 rtr =
glsc3(tmp, coef%mult, coef%binv, n)
259 rnorm = sqrt(rtr*norm_fac)
260 ksp_results%res_start = rnorm
261 ksp_results%res_final = rnorm
263 if (rnorm .eq. zero)
return
265 call this%monitor_start(
'cpldCG')
266 do iter = 1, max_iter
267 call this%M%solve(z1, this%r1, n)
268 call this%M%solve(z2, this%r2, n)
269 call this%M%solve(z3, this%r3, n)
272 do concurrent(i = 1:n)
273 this%tmp(i) = z1(i) * r1(i) &
278 rtz1 =
glsc2(tmp, coef%mult, n)
281 if (iter .eq. 1) beta = zero
282 do concurrent(i = 1:n)
283 p1(i) = p1(i) * beta + z1(i)
284 p2(i) = p2(i) * beta + z2(i)
285 p3(i) = p3(i) * beta + z3(i)
288 call ax%compute_vector(w1, w2, w3, p1, p2, p3, coef, x%msh, x%Xh)
289 call gs_h%op(w1, n, gs_op_add)
290 call gs_h%op(w2, n, gs_op_add)
291 call gs_h%op(w3, n, gs_op_add)
292 call blstx%apply_scalar(w1, n)
293 call blsty%apply_scalar(w2, n)
294 call blstz%apply_scalar(w3, n)
296 do concurrent(i = 1:n)
297 tmp(i) = w1(i) * p1(i) &
302 pap =
glsc2(tmp, coef%mult, n)
306 do concurrent(i = 1:n)
307 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * p1(i)
308 y%x(i,1,1,1) = y%x(i,1,1,1) + alpha * p2(i)
309 z%x(i,1,1,1) = z%x(i,1,1,1) + alpha * p3(i)
310 r1(i) = r1(i) + alphm * w1(i)
311 r2(i) = r2(i) + alphm * w2(i)
312 r3(i) = r3(i) + alphm * w3(i)
313 tmp(i) = r1(i)**2 + r2(i)**2 + r3(i)**2
316 rtr =
glsc3(tmp, coef%mult, coef%binv, n)
317 if (iter .eq. 1) rtr0 = rtr
318 rnorm = sqrt(rtr * norm_fac)
319 call this%monitor_iter(iter, rnorm)
320 if (rnorm .lt. this%abs_tol)
then
325 call this%monitor_stop()
326 ksp_results%res_final = rnorm
327 ksp_results%iter = iter
328 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.