162 function pipecg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter)
result(ksp_results)
163 class(
pipecg_t),
intent(inout) :: this
164 class(
ax_t),
intent(in) :: ax
165 type(
field_t),
intent(inout) :: x
166 integer,
intent(in) :: n
167 real(kind=
rp),
dimension(n),
intent(in) :: f
168 type(
coef_t),
intent(inout) :: coef
170 type(
gs_t),
intent(inout) :: gs_h
172 integer,
optional,
intent(in) :: niter
173 integer :: iter, max_iter, i, j, k, ierr, p_cur, p_prev, u_prev
174 real(kind=
rp) :: rnorm, rtr, reduction(3), norm_fac
176 real(kind=
rp) :: gamma1, gamma2, delta
178 type(mpi_request) :: request
179 type(mpi_status) :: status
181 if (
present(niter))
then
184 max_iter = this%max_iter
186 norm_fac = 1.0_rp / sqrt(coef%volume)
188 associate(p => this%p, q => this%q, r => this%r, s => this%s, &
189 u => this%u, w => this%w, z => this%z, mi => this%mi, ni => this%ni)
200 call this%M%solve(u(1,u_prev), r, n)
201 call ax%compute(w, u(1,u_prev), coef, x%msh, x%Xh)
202 call gs_h%op(w, n, gs_op_add)
203 call blst%apply(w, n)
205 rtr =
glsc3(r, coef%mult, r, n)
206 rnorm = sqrt(rtr)*norm_fac
207 ksp_results%res_start = rnorm
208 ksp_results%res_final = rnorm
211 if(
abscmp(rnorm, 0.0_rp))
then
212 ksp_results%converged = .true.
221 tmp1 = tmp1 + r(i) * coef%mult(i,1,1,1) * u(i,u_prev)
222 tmp2 = tmp2 + w(i) * coef%mult(i,1,1,1) * u(i,u_prev)
223 tmp3 = tmp3 + r(i) * coef%mult(i,1,1,1) * r(i)
229 call this%monitor_start(
'PipeCG')
230 do iter = 1, max_iter
231 call mpi_iallreduce(mpi_in_place, reduction, 3, &
234 call this%M%solve(mi, w, n)
235 call ax%compute(ni, mi, coef, x%msh, x%Xh)
236 call gs_h%op(ni, n, gs_op_add)
237 call blst%apply(ni, n)
239 call mpi_wait(request, status, ierr)
241 gamma1 = reduction(1)
245 rnorm = sqrt(rtr)*norm_fac
246 call this%monitor_iter(iter, rnorm)
247 if (rnorm .lt. this%abs_tol)
exit
249 if (iter .gt. 1)
then
250 beta(p_cur) = gamma1 / gamma2
251 alpha(p_cur) = gamma1 / (delta - (beta(p_cur) * gamma1/alpha(p_prev)))
254 alpha(p_cur) = gamma1/delta
263 z(i+k) = beta(p_cur) * z(i+k) + ni(i+k)
264 q(i+k) = beta(p_cur) * q(i+k) + mi(i+k)
265 s(i+k) = beta(p_cur) * s(i+k) + w(i+k)
266 r(i+k) = r(i+k) - alpha(p_cur) * s(i+k)
267 u(i+k,p_cur) = u(i+k,u_prev) - alpha(p_cur) * q(i+k)
268 w(i+k) = w(i+k) - alpha(p_cur) * z(i+k)
269 tmp1 = tmp1 + r(i+k) * coef%mult(i+k,1,1,1) * u(i+k,p_cur)
270 tmp2 = tmp2 + w(i+k) * coef%mult(i+k,1,1,1) * u(i+k,p_cur)
271 tmp3 = tmp3 + r(i+k) * coef%mult(i+k,1,1,1) * r(i+k)
275 z(i+k) = beta(p_cur) * z(i+k) + ni(i+k)
276 q(i+k) = beta(p_cur) * q(i+k) + mi(i+k)
277 s(i+k) = beta(p_cur) * s(i+k) + w(i+k)
278 r(i+k) = r(i+k) - alpha(p_cur) * s(i+k)
279 u(i+k,p_cur) = u(i+k,u_prev) - alpha(p_cur) * q(i+k)
280 w(i+k) = w(i+k) - alpha(p_cur) * z(i+k)
281 tmp1 = tmp1 + r(i+k) * coef%mult(i+k,1,1,1) * u(i+k,p_cur)
282 tmp2 = tmp2 + w(i+k) * coef%mult(i+k,1,1,1) * u(i+k,p_cur)
283 tmp3 = tmp3 + r(i+k) * coef%mult(i+k,1,1,1) * r(i+k)
301 p(i+k) = beta(j) * p(i+k) + u(i+k,p_prev)
302 x_plus(k) = x_plus(k) + alpha(j) * p(i+k)
307 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(k)
315 p(i+k) = beta(j) * p(i+k) + u(i+k,p_prev)
316 x_plus(1) = x_plus(1) + alpha(j) * p(i+k)
319 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(1)
326 alpha(1) = alpha(p_cur)
327 beta(1) = beta(p_cur)
336 if ( p_cur .ne. 1)
then
345 p(i+k) = beta(j) * p(i+k) + u(i+k,p_prev)
346 x_plus(k) = x_plus(k) + alpha(j) * p(i+k)
351 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(k)
359 p(i+k) = beta(j) * p(i+k) + u(i+k,p_prev)
360 x_plus(1) = x_plus(1) + alpha(j) * p(i+k)
363 x%x(i+k,1,1,1) = x%x(i+k,1,1,1) + x_plus(1)
369 call this%monitor_stop()
370 ksp_results%res_final = rnorm
371 ksp_results%iter = iter
372 ksp_results%converged = this%is_converged(iter, rnorm)
380 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
381 class(
pipecg_t),
intent(inout) :: this
382 class(ax_t),
intent(in) :: ax
383 type(field_t),
intent(inout) :: x
384 type(field_t),
intent(inout) :: y
385 type(field_t),
intent(inout) :: z
386 integer,
intent(in) :: n
387 real(kind=rp),
dimension(n),
intent(in) :: fx
388 real(kind=rp),
dimension(n),
intent(in) :: fy
389 real(kind=rp),
dimension(n),
intent(in) :: fz
390 type(coef_t),
intent(inout) :: coef
391 type(bc_list_t),
intent(inout) :: blstx
392 type(bc_list_t),
intent(inout) :: blsty
393 type(bc_list_t),
intent(inout) :: blstz
394 type(gs_t),
intent(inout) :: gs_h
395 type(ksp_monitor_t),
dimension(3) :: ksp_results
396 integer,
optional,
intent(in) :: niter
398 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
399 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
400 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)