54 real(kind=
rp),
allocatable :: w(:)
55 real(kind=
rp),
allocatable :: r(:)
56 real(kind=
rp),
allocatable :: z(:)
57 real(kind=
rp),
allocatable :: p(:,:)
58 real(kind=
rp),
allocatable :: alpha(:)
59 type(c_ptr) :: w_d = c_null_ptr
60 type(c_ptr) :: r_d = c_null_ptr
61 type(c_ptr) :: z_d = c_null_ptr
62 type(c_ptr) :: alpha_d = c_null_ptr
63 type(c_ptr) :: p_d_d = c_null_ptr
64 type(c_ptr),
allocatable :: p_d(:)
65 type(c_ptr) :: gs_event = c_null_ptr
76 bind(c, name=
'cuda_fusedcg_update_p')
77 use,
intrinsic :: iso_c_binding
80 type(c_ptr),
value :: p_d, z_d, po_d
88 bind(c, name=
'cuda_fusedcg_update_x')
89 use,
intrinsic :: iso_c_binding
91 type(c_ptr),
value :: x_d, p_d, alpha
92 integer(c_int) :: p_cur, n
98 p_cur, n) bind(c, name=
'cuda_fusedcg_part2')
99 use,
intrinsic :: iso_c_binding
102 type(c_ptr),
value :: a_d, b_d, c_d, alpha_d
104 integer(c_int) :: n, p_cur
110 bind(c, name=
'hip_fusedcg_update_p')
111 use,
intrinsic :: iso_c_binding
114 type(c_ptr),
value :: p_d, z_d, po_d
122 bind(c, name=
'hip_fusedcg_update_x')
123 use,
intrinsic :: iso_c_binding
125 type(c_ptr),
value :: x_d, p_d, alpha
126 integer(c_int) :: p_cur, n
132 p_cur, n) bind(c, name=
'hip_fusedcg_part2')
133 use,
intrinsic :: iso_c_binding
136 type(c_ptr),
value :: a_d, b_d, c_d, alpha_d
138 integer(c_int) :: n, p_cur
146 type(c_ptr),
value :: p_d, z_d, po_d
154 call neko_error(
'No device backend configured')
159 type(c_ptr),
value :: x_d, p_d, alpha
160 integer(c_int) :: p_cur, n
166 call neko_error(
'No device backend configured')
171 p_cur, n)
result(res)
172 type(c_ptr),
value :: a_d, b_d, c_d, alpha_d
182 call neko_error(
'No device backend configured')
185#ifndef HAVE_DEVICE_MPI
186 if (pe_size .gt. 1)
then
187 call mpi_allreduce(mpi_in_place, res, 1, &
188 mpi_real_precision, mpi_sum, neko_comm, ierr)
198 class(pc_t),
optional,
intent(in),
target :: M
199 integer,
intent(in) :: n
200 integer,
intent(in) :: max_iter
201 real(kind=rp),
optional,
intent(in) :: rel_tol
202 real(kind=rp),
optional,
intent(in) :: abs_tol
203 logical,
optional,
intent(in) :: monitor
205 integer(c_size_t) :: p_size
221 call device_map(this%w, this%w_d, n)
222 call device_map(this%r, this%r_d, n)
223 call device_map(this%z, this%z_d, n)
226 this%p_d(i) = c_null_ptr
227 call device_map(this%p(:,i), this%p_d(i), n)
231 call device_alloc(this%p_d_d, p_size)
232 ptr = c_loc(this%p_d)
233 call device_memcpy(ptr, this%p_d_d, p_size, &
234 host_to_device, sync=.false.)
235 if (
present(rel_tol) .and.
present(abs_tol) .and.
present(monitor))
then
236 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
237 else if (
present(rel_tol) .and.
present(abs_tol))
then
238 call this%ksp_init(max_iter, rel_tol, abs_tol)
239 else if (
present(monitor) .and.
present(abs_tol))
then
240 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
241 else if (
present(rel_tol) .and.
present(monitor))
then
242 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
243 else if (
present(rel_tol))
then
244 call this%ksp_init(max_iter, rel_tol = rel_tol)
245 else if (
present(abs_tol))
then
246 call this%ksp_init(max_iter, abs_tol = abs_tol)
247 else if (
present(monitor))
then
248 call this%ksp_init(max_iter, monitor = monitor)
250 call this%ksp_init(max_iter)
253 call device_event_create(this%gs_event, 2)
264 if (
allocated(this%w))
then
268 if (
allocated(this%r))
then
273 if (
allocated(this%z))
then
278 if (
allocated(this%alpha))
then
279 deallocate(this%alpha)
282 if (
allocated(this%p))
then
286 if (c_associated(this%w_d))
then
287 call device_free(this%w_d)
290 if (c_associated(this%r_d))
then
291 call device_free(this%r_d)
294 if (c_associated(this%z_d))
then
295 call device_free(this%z_d)
298 if (c_associated(this%alpha_d))
then
299 call device_free(this%alpha_d)
302 if (
allocated(this%p_d))
then
304 if (c_associated(this%p_d(i)))
then
305 call device_free(this%p_d(i))
312 if (c_associated(this%gs_event))
then
313 call device_event_destroy(this%gs_event)
321 class(ax_t),
intent(in) :: ax
322 type(field_t),
intent(inout) :: x
323 integer,
intent(in) :: n
324 real(kind=rp),
dimension(n),
intent(in) :: f
325 type(coef_t),
intent(inout) :: coef
326 type(bc_list_t),
intent(inout) :: blst
327 type(gs_t),
intent(inout) :: gs_h
328 type(ksp_monitor_t) :: ksp_results
329 integer,
optional,
intent(in) :: niter
330 integer :: iter, max_iter, ierr, i, p_cur, p_prev
331 real(kind=rp) :: rnorm, rtr, norm_fac, rtz1, rtz2
332 real(kind=rp) :: pap, beta
334 f_d = device_get_ptr(f)
336 if (
present(niter))
then
339 max_iter = ksp_max_iter
341 norm_fac = 1.0_rp / sqrt(coef%volume)
343 associate(w => this%w, r => this%r, p => this%p, z => this%z, &
344 alpha => this%alpha, alpha_d => this%alpha_d, &
345 w_d => this%w_d, r_d => this%r_d, z_d => this%z_d, &
346 p_d => this%p_d, p_d_d => this%p_d_d)
351 call device_rzero(x%x_d, n)
352 call device_rzero(p_d(1), n)
353 call device_copy(r_d, f_d, n)
355 rtr = device_glsc3(r_d, coef%mult_d, r_d, n)
356 rnorm = sqrt(rtr)*norm_fac
357 ksp_results%res_start = rnorm
358 ksp_results%res_final = rnorm
360 if(abscmp(rnorm, 0.0_rp))
return
362 call this%monitor_start(
'FusedCG')
363 do iter = 1, max_iter
364 call this%M%solve(z, r, n)
366 rtz1 = device_glsc3(r_d, coef%mult_d, z_d, n)
368 if (iter .eq. 1) beta = 0.0_rp
371 call ax%compute(w, p(1, p_cur), coef, x%msh, x%Xh)
372 call gs_h%op(w, n, gs_op_add, this%gs_event)
373 call device_event_sync(this%gs_event)
374 call blst%apply(w, n)
376 pap = device_glsc3(w_d, coef%mult_d, this%p_d(p_cur), n)
378 alpha(p_cur) = rtz1 / pap
380 alpha_d, alpha(p_cur), p_cur, n)
381 rnorm = sqrt(rtr)*norm_fac
382 call this%monitor_iter(iter, rnorm)
384 (rnorm .lt. this%abs_tol) .or. iter .eq. max_iter)
then
388 if (rnorm .lt. this%abs_tol)
exit
394 call this%monitor_stop()
395 ksp_results%res_final = rnorm
396 ksp_results%iter = iter
397 ksp_results%converged = this%is_converged(iter, rnorm)
405 n, coef, blstx, blsty, blstz, gs_h, niter)
result(ksp_results)
407 class(ax_t),
intent(in) :: ax
408 type(field_t),
intent(inout) :: x
409 type(field_t),
intent(inout) :: y
410 type(field_t),
intent(inout) :: z
411 integer,
intent(in) :: n
412 real(kind=rp),
dimension(n),
intent(in) :: fx
413 real(kind=rp),
dimension(n),
intent(in) :: fy
414 real(kind=rp),
dimension(n),
intent(in) :: fz
415 type(coef_t),
intent(inout) :: coef
416 type(bc_list_t),
intent(inout) :: blstx
417 type(bc_list_t),
intent(inout) :: blsty
418 type(bc_list_t),
intent(inout) :: blstz
419 type(gs_t),
intent(inout) :: gs_h
420 type(ksp_monitor_t),
dimension(3) :: ksp_results
421 integer,
optional,
intent(in) :: niter
423 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
424 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
425 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
void hip_fusedcg_update_x(void *x, void *p, void *alpha, int *p_cur, int *n)
real hip_fusedcg_part2(void *a, void *b, void *c, void *alpha_d, real *alpha, int *p_cur, int *n)
void hip_fusedcg_update_p(void *p, void *z, void *po, real *beta, int *n)
Defines a Matrix-vector product.
subroutine, public device_rzero(a_d, n)
Zero a real vector.
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n)
Weighted inner product .
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
Device abstraction, common interface for various accelerators.
Defines a fused Conjugate Gradient method for accelerators.
subroutine fusedcg_device_free(this)
Deallocate a pipelined PCG solver.
real(kind=rp) function device_fusedcg_part2(a_d, b_d, c_d, alpha_d, alpha, p_cur, n)
type(ksp_monitor_t) function, dimension(3) fusedcg_device_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Pipelined PCG solve coupled solve.
integer, parameter device_fusedcg_p_space
subroutine device_fusedcg_update_x(x_d, p_d, alpha, p_cur, n)
subroutine fusedcg_device_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a fused PCG solver.
type(ksp_monitor_t) function fusedcg_device_solve(this, ax, x, f, n, coef, blst, gs_h, niter)
Pipelined PCG solve.
subroutine device_fusedcg_update_p(p_d, z_d, po_d, beta, n)
Implements the base abstract type for Krylov solvers plus helper types.
integer, parameter, public ksp_max_iter
Maximum number of iters.
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.
Base type for a matrix-vector product providing .
A list of allocatable `bc_t`. Follows the standard interface of lists.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Fused preconditioned conjugate gradient method.
Type for storing initial and final residuals in a Krylov solver.
Base abstract type for a canonical Krylov method, solving .
Defines a canonical Krylov preconditioner.