85 use,
intrinsic :: iso_c_binding, only : c_ptr, c_size_t, &
86 c_sizeof, c_null_ptr, c_loc, c_associated
92 real(kind=
rp),
allocatable :: xx(:,:)
93 real(kind=
rp),
allocatable :: bb(:,:)
94 real(kind=
rp),
allocatable :: xbar(:)
95 type(c_ptr),
allocatable :: xx_d(:)
96 type(c_ptr),
allocatable :: bb_d(:)
97 type(c_ptr) :: xbar_d = c_null_ptr
98 type(c_ptr) :: alpha_d = c_null_ptr
99 type(c_ptr) :: xx_d_d = c_null_ptr
100 type(c_ptr) :: bb_d_d = c_null_ptr
102 real(kind=
rp) :: tol = 1e-7_rp
104 real(kind=
rp) :: proj_res
105 integer :: proj_m = 0
106 integer :: activ_step
122 integer,
intent(in) :: n
123 integer,
optional,
intent(in) :: L, activ_step
125 integer(c_size_t) :: ptr_size
137 if (
present(activ_step))
then
138 this%activ_step = activ_step
145 allocate(this%xx(n, this%L))
146 allocate(this%bb(n, this%L))
147 allocate(this%xbar(n))
148 allocate(this%xx_d(this%L))
149 allocate(this%bb_d(this%L))
150 call rzero(this%xbar, n)
152 call rzero(this%xx(1, i), n)
153 call rzero(this%bb(1, i), n)
158 call device_alloc(this%alpha_d, int(c_sizeof(dummy)*this%L, c_size_t))
164 this%xx_d(i) = c_null_ptr
165 call device_map(this%xx(:, i), this%xx_d(i), n)
167 this%bb_d(i) = c_null_ptr
168 call device_map(this%bb(:, i), this%bb_d(i), n)
172 ptr_size = c_sizeof(c_null_ptr) * this%L
174 ptr = c_loc(this%xx_d)
178 ptr = c_loc(this%bb_d)
189 if (
allocated(this%xx))
then
192 if (
allocated(this%bb))
then
195 if (
allocated(this%xbar))
then
196 deallocate(this%xbar)
198 if (
allocated(this%xx_d))
then
200 if (c_associated(this%xx_d(i)))
then
205 if (c_associated(this%xx_d_d))
then
208 if (c_associated(this%xbar_d))
then
211 if (c_associated(this%alpha_d))
then
214 if (
allocated(this%bb_d))
then
216 if (c_associated(this%bb_d(i)))
then
221 if (c_associated(this%bb_d_d))
then
230 integer,
intent(inout) :: n
231 real(kind=
rp),
intent(inout),
dimension(n) :: b
232 integer,
intent(in) :: tstep
233 class(
coef_t),
intent(inout) :: coef
235 character(len=*),
optional :: string
237 if (tstep .gt. this%activ_step .and. this%L .gt. 0)
then
238 if (dt_controller%if_variable_dt)
then
240 if (dt_controller%dt_last_change .eq. 0)
then
242 else if (dt_controller%dt_last_change .gt. this%activ_step - 1)
then
245 call this%project_on(b, coef, n)
246 if (
present(string))
then
247 call this%log_info(string)
251 call this%project_on(b, coef, n)
252 if (
present(string))
then
253 call this%log_info(string)
263 integer,
intent(inout) :: n
264 class(
ax_t),
intent(inout) :: Ax
265 class(
coef_t),
intent(inout) :: coef
267 type(
gs_t),
intent(inout) :: gs_h
268 real(kind=
rp),
intent(inout),
dimension(n) :: x
269 integer,
intent(in) :: tstep
272 if (tstep .gt. this%activ_step .and. this%L .gt. 0)
then
273 if (.not.(dt_controller%if_variable_dt) .or. &
274 (dt_controller%dt_last_change .gt. this%activ_step - 1))
then
275 call this%project_back(x, ax, coef, bclst, gs_h, n)
283 integer,
intent(inout) :: n
284 class(
coef_t),
intent(inout) :: coef
285 real(kind=
rp),
intent(inout),
dimension(n) :: b
297 integer,
intent(inout) :: n
298 class(
ax_t),
intent(inout) :: Ax
299 class(
coef_t),
intent(inout) :: coef
301 type(
gs_t),
intent(inout) :: gs_h
302 real(kind=
rp),
intent(inout),
dimension(n) :: x
310 if (this%m .gt. 0)
call device_add2(x_d, this%xbar_d, n)
311 if (this%m .eq. this%L)
then
314 this%m = min(this%m+1, this%L)
320 if (this%m .gt. 0)
call add2(x, this%xbar, n)
321 if (this%m .eq. this%L)
then
324 this%m = min(this%m+1, this%L)
327 call copy(this%xx(1, this%m), x, n)
330 call ax%compute(this%bb(1, this%m), x, coef, coef%msh, coef%Xh)
331 call gs_h%gs_op_vector(this%bb(1, this%m), n, gs_op_add)
332 call bclst%apply_scalar(this%bb(1, this%m), n)
342 integer,
intent(inout) :: n
343 class(
coef_t),
intent(inout) :: coef
344 real(kind=
rp),
intent(inout),
dimension(n) :: b
345 integer :: i, j, k, l, ierr
346 real(kind=
rp) :: work(this%L), alpha(this%L), s
348 associate(xbar => this%xbar, xx => this%xx, &
351 if (this%m .le. 0)
return
354 call rzero(alpha, this%m)
355 this%proj_res = sqrt(
glsc3(b, b, coef%mult, n) / coef%volume)
362 s = s + xx(i+l, k) * coef%mult(i+l,1,1,1) * b(i+l)
364 alpha(k) = alpha(k) + s
369 call mpi_allreduce(mpi_in_place, alpha, this%m, &
372 call rzero(work, this%m)
377 xbar(i+l) = alpha(1) * xx(i+l,1)
378 b(i+l) = b(i+l) - alpha(1) * bb(i+l,1)
382 xbar(i+l) = xbar(i+l) + alpha(k) * xx(i+l,k)
383 b(i+l) = b(i+l)- alpha(k) * bb(i+l,k)
390 s = s + xx(i+l,k) * coef%mult(i+l,1,1,1) * b(i+l)
392 work(k) = work(k) + s
396 call mpi_allreduce(work, alpha, this%m, &
403 xbar(i+l) = xbar(i+l) + alpha(k) * xx(i+l,k)
404 b(i+l) = b(i+l) - alpha(k) * bb(i+l,k)
413 integer,
intent(inout) :: n
414 class(coef_t),
intent(inout) :: coef
415 real(kind=rp),
intent(inout),
dimension(n) :: b
416 real(kind=rp) :: alpha(this%L)
419 b_d = device_get_ptr(b)
421 associate(xbar_d => this%xbar_d, xx_d => this%xx_d, xx_d_d => this%xx_d_d, &
422 bb_d => this%bb_d, bb_d_d => this%bb_d_d, alpha_d => this%alpha_d)
424 if (this%m .le. 0)
return
428 this%proj_res = sqrt(device_glsc3(b_d, b_d, coef%mult_d, n)/coef%volume)
430 if (neko_device_mpi .and. (neko_bcknd_opencl .ne. 1))
then
431 call device_proj_on(alpha_d, b_d, xx_d_d, bb_d_d, &
432 coef%mult_d, xbar_d, this%m, n)
434 if (neko_bcknd_opencl .eq. 1)
then
436 alpha(i) = device_glsc3(b_d, xx_d(i), coef%mult_d, n)
439 call device_glsc3_many(alpha, b_d, xx_d_d, coef%mult_d, this%m, n)
440 call device_memcpy(alpha, alpha_d, this%m, &
441 host_to_device, sync = .false.)
443 call device_rzero(xbar_d, n)
444 if (neko_bcknd_opencl .eq. 1)
then
446 call device_add2s2(xbar_d, xx_d(i), alpha(i), n)
448 call cmult(alpha, -1.0_rp, this%m)
450 call device_add2s2_many(xbar_d, xx_d_d, alpha_d, this%m, n)
451 call device_cmult(alpha_d, -1.0_rp, this%m)
454 if (neko_bcknd_opencl .eq. 1)
then
456 call device_add2s2(b_d, bb_d(i), alpha(i), n)
457 alpha(i) = device_glsc3(b_d, xx_d(i), coef%mult_d, n)
460 call device_add2s2_many(b_d, bb_d_d, alpha_d, this%m, n)
461 call device_glsc3_many(alpha, b_d, xx_d_d, coef%mult_d, this%m, n)
462 call device_memcpy(alpha, alpha_d, this%m, &
463 host_to_device, sync = .false.)
466 if (neko_bcknd_opencl .eq. 1)
then
468 call device_add2s2(xbar_d, xx_d(i), alpha(i), n)
469 call cmult(alpha, -1.0_rp, this%m)
470 call device_add2s2(b_d, bb_d(i), alpha(i), n)
473 call device_add2s2_many(xbar_d, xx_d_d, alpha_d, this%m, n)
474 call device_cmult(alpha_d, -1.0_rp, this%m)
475 call device_add2s2_many(b_d, bb_d_d, alpha_d, this%m, n)
485 type(coef_t),
intent(in) :: coef
486 integer,
intent(in) :: n
488 if (neko_bcknd_device .eq. 1)
then
498 integer,
intent(in) :: n
499 type(c_ptr),
dimension(this%L) :: xx_d, bb_d
500 type(c_ptr),
intent(in) :: w_d
501 real(kind=rp) :: nrm, scl
502 real(kind=rp) :: alpha(this%L)
505 associate(m => this%m, xx_d_d => this%xx_d_d, &
506 bb_d_d => this%bb_d_d, alpha_d => this%alpha_d)
510 if (neko_device_mpi .and. (neko_bcknd_opencl .ne. 1))
then
511 call device_project_ortho(alpha_d, bb_d(m), xx_d_d, bb_d_d, &
512 w_d, xx_d(m), this%m, n, nrm)
514 if (neko_bcknd_opencl .eq. 1)
then
516 alpha(i) = device_glsc3(bb_d(m), xx_d(i), w_d,n)
519 call device_glsc3_many(alpha, bb_d(m), xx_d_d, w_d, m, n)
522 call cmult(alpha, -1.0_rp,m)
523 if (neko_bcknd_opencl .eq. 1)
then
525 call device_add2s2(xx_d(m), xx_d(i), alpha(i), n)
526 call device_add2s2(bb_d(m), bb_d(i), alpha(i), n)
528 alpha(i) = device_glsc3(bb_d(m), xx_d(i), w_d, n)
531 call device_memcpy(alpha, alpha_d, this%m, &
532 host_to_device, sync = .false.)
533 call device_add2s2_many(xx_d(m), xx_d_d, alpha_d, m-1, n)
534 call device_add2s2_many(bb_d(m), bb_d_d, alpha_d, m-1, n)
536 call device_glsc3_many(alpha, bb_d(m), xx_d_d, w_d, m, n)
538 call cmult(alpha, -1.0_rp,m)
539 if (neko_bcknd_opencl .eq. 1)
then
541 call device_add2s2(xx_d(m), xx_d(i), alpha(i), n)
542 call device_add2s2(bb_d(m), bb_d(i), alpha(i), n)
543 alpha(i) = device_glsc3(bb_d(m), xx_d(i), w_d, n)
546 call device_memcpy(alpha, alpha_d, m, &
547 host_to_device, sync = .false.)
548 call device_add2s2_many(xx_d(m), xx_d_d, alpha_d, m-1, n)
549 call device_add2s2_many(bb_d(m), bb_d_d, alpha_d, m-1, n)
550 call device_glsc3_many(alpha, bb_d(m), xx_d_d, w_d, m, n)
554 alpha(m) = device_glsc3(xx_d(m), w_d, bb_d(m), n)
555 alpha(m) = sqrt(alpha(m))
557 if (alpha(m) .gt. this%tol*nrm)
then
558 scl = 1.0_rp / alpha(m)
559 call device_cmult(xx_d(m), scl, n)
560 call device_cmult(bb_d(m), scl, n)
564 if (pe_rank .eq. 0)
then
565 call neko_warning(
'New vector not linearly indepependent!')
577 integer,
intent(in) :: n
578 real(kind=rp),
dimension(n, this%L),
intent(inout) :: xx, bb
579 real(kind=rp),
dimension(n),
intent(in) :: w
580 real(kind=rp) :: nrm, scl1, scl2, c, s
581 real(kind=rp) :: alpha(this%L), beta(this%L)
582 integer :: i, j, k, l, h, ierr
584 associate(m => this%m)
591 do i = 1, n, neko_blk_size
592 j = min(neko_blk_size, n-i+1)
597 s = s + xx(i+l,k) * w(i+l) * bb(i+l,m)
598 c = c + bb(i+l,k) * w(i+l) * xx(i+l,m)
600 alpha(k) = alpha(k) + 0.5_rp * (s + c)
604 call mpi_allreduce(mpi_in_place, alpha, this%m, &
605 mpi_real_precision, mpi_sum, neko_comm, ierr)
610 do i = 1, n, neko_blk_size
611 j = min(neko_blk_size, n-i+1)
614 xx(i+l,m) = xx(i+l,m) - alpha(k) * xx(i+l,k)
615 bb(i+l,m) = bb(i+l,m) - alpha(k) * bb(i+l,k)
621 do i = 1, n, neko_blk_size
622 j = min(neko_blk_size, n-i+1)
627 s = s + xx(i+l,k) * w(i+l) * bb(i+l,m)
628 c = c + bb(i+l,k) * w(i+l) * xx(i+l,m)
630 beta(k) = beta(k) + 0.5_rp * (s + c)
634 call mpi_allreduce(mpi_in_place, beta, this%m-1, &
635 mpi_real_precision, mpi_sum, neko_comm, ierr)
639 do i = 1, n, neko_blk_size
640 j = min(neko_blk_size,n-i+1)
643 xx(i+l,m) = xx(i+l,m) - beta(k) * xx(i+l,k)
644 bb(i+l,m) = bb(i+l,m) - beta(k) * bb(i+l,k)
649 s = s + xx(i+l,m) * w(i+l) * bb(i+l,m)
651 alpha(m) = alpha(m) + s
654 alpha(k) = alpha(k) + beta(k)
658 call mpi_allreduce(mpi_in_place, alpha(m), 1, &
659 mpi_real_precision, mpi_sum, neko_comm, ierr)
660 alpha(m) = sqrt(alpha(m))
663 if (alpha(m) .gt. this%tol*nrm)
then
665 scl1 = 1.0_rp / alpha(m)
667 xx(1+i,m) = scl1 * xx(1+i,m)
668 bb(1+i,m) = scl1 * bb(1+i,m)
673 if (pe_rank .eq. 0)
then
674 call neko_warning(
'New vector not linearly indepependent!')
685 character(len=*),
intent(in) :: string
686 character(len=LOG_SIZE) :: log_buf
688 if (this%proj_m .gt. 0)
then
689 write(log_buf,
'(A,A)')
'Projection ', string
690 call neko_log%message(log_buf)
691 write(log_buf,
'(A,A)')
'Proj. vec.:',
' Orig. residual:'
692 call neko_log%message(log_buf)
693 write(log_buf,
'(I11,3x, E15.7,5x)') this%proj_m, this%proj_res
694 call neko_log%message(log_buf)
701 integer,
intent(in) :: n
708 if (neko_bcknd_device .eq. 1)
then
709 call device_rzero(this%xx_d(i), n)
710 call device_rzero(this%xx_d(i), n)
713 this%xx(j,i) = 0.0_rp
714 this%bb(j,i) = 0.0_rp
Return the device pointer for an associated Fortran array.
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Defines a Matrix-vector product.
type(mpi_comm) neko_comm
MPI communicator.
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
subroutine, public device_add2(a_d, b_d, n)
Vector addition .
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_add2s2(a_d, b_d, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
subroutine, public device_add2s2_many(y_d, x_d_d, a_d, j, n)
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 .
subroutine, public device_glsc3_many(h, w_d, v_d_d, mult_d, j, n)
Interface for device projection.
subroutine, public device_proj_on(alpha_d, b_d, x_d_d, b_d_d, mult_d, xbar_d, j, n)
subroutine, public device_project_ortho(alpha_d, b_d, x_d_d, b_d_d, w_d, xm_d, j, n, nrm)
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public cmult(a, c, n)
Multiplication by constant c .
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
subroutine, public add2(a, b, n)
Vector addition .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_blk_size
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_opencl
logical, parameter neko_device_mpi
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Project x onto X, the space of old solutions and back again.
subroutine device_project_on(this, b, coef, n)
subroutine projection_free(this)
subroutine bcknd_clear(this, n)
subroutine device_proj_ortho(this, xx_d, bb_d, w_d, n)
subroutine bcknd_project_back(this, x, ax, coef, bclst, gs_h, n)
subroutine cpu_proj_ortho(this, xx, bb, w, n)
subroutine, public proj_ortho(this, coef, n)
subroutine bcknd_project_on(this, b, coef, n)
subroutine projection_init(this, n, l, activ_step)
subroutine cpu_project_on(this, b, coef, n)
subroutine projection_pre_solving(this, b, tstep, coef, n, dt_controller, string)
subroutine projection_post_solving(this, x, ax, coef, bclst, gs_h, n, tstep, dt_controller)
subroutine print_proj_info(this, string)
Implements type time_step_controller.
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
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,...
Provides a tool to set time step dt.