44  use, 
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
 
   50     real(kind=
rp), 
allocatable :: d(:,:,:,:)
 
   51     type(c_ptr) :: d_d = c_null_ptr
 
   52     type(
gs_t), 
pointer :: gs_h
 
   55     type(c_ptr) :: gs_event = c_null_ptr
 
 
   65          G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
 
   66          bind(c, name=
'hip_jacobi_update')
 
   67       use, 
intrinsic :: iso_c_binding
 
   68       type(c_ptr), 
value :: d_d, dxt_d, dyt_d, dzt_d
 
   69       type(c_ptr), 
value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
 
   70       integer(c_int) :: nelv, lx
 
 
   76          G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
 
   77          bind(c, name=
'cuda_jacobi_update')
 
   78       use, 
intrinsic :: iso_c_binding
 
   79       type(c_ptr), 
value :: d_d, dxt_d, dyt_d, dzt_d
 
   80       type(c_ptr), 
value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
 
   81       integer(c_int) :: nelv, lx
 
 
   87          G11_d, G22_d, G33_d, G12_d, G13_d, G23_d, nelv, lx) &
 
   88          bind(c, name=
'opencl_jacobi_update')
 
   89       use, 
intrinsic :: iso_c_binding
 
   90       type(c_ptr), 
value :: d_d, dxt_d, dyt_d, dzt_d
 
   91       type(c_ptr), 
value :: G11_d, G22_d, G33_d, G12_d, G13_d, G23_d
 
   92       integer(c_int) :: nelv, lx
 
 
  100    type(
coef_t), 
intent(in), 
target :: coef
 
  101    type(
dofmap_t), 
intent(in), 
target :: dof
 
  102    type(
gs_t), 
intent(inout), 
target :: gs_h
 
  110    allocate(this%d(dof%Xh%lx,dof%Xh%ly,dof%Xh%lz, dof%msh%nelv))
 
  112    call device_map(this%d, this%d_d, 
size(this%d))
 
 
  123    if (c_associated(this%d_d)) 
then 
  125       this%d_d = c_null_ptr
 
  128    if (
allocated(this%d)) 
then 
  132    if (c_associated(this%gs_event)) 
then 
 
  144    integer, 
intent(in) :: n
 
  146    real(kind=
rp), 
dimension(n), 
intent(inout) :: z
 
  147    real(kind=
rp), 
dimension(n), 
intent(inout) :: r
 
  148    type(c_ptr) :: z_d, r_d
 
 
  159    integer :: lz, ly, lx
 
  160    associate(dof => this%dof, coef => this%coef, xh => this%dof%Xh, &
 
  161         gs_h => this%gs_h, nelv => this%dof%msh%nelv)
 
  169           coef%G11_d, coef%G22_d, coef%G33_d, &
 
  170           coef%G12_d, coef%G13_d, coef%G23_d, &
 
  174           coef%G11_d, coef%G22_d, coef%G33_d, &
 
  175           coef%G12_d, coef%G13_d, coef%G23_d, &
 
  179           coef%G11_d, coef%G22_d, coef%G33_d, &
 
  180           coef%G12_d, coef%G13_d, coef%G23_d, &
 
  184      call device_col2(this%d_d, coef%h1_d, coef%dof%size())
 
  187         call device_addcol3(this%d_d, coef%h2_d, coef%B_d, coef%dof%size())
 
  190      call gs_h%op(this%d, dof%size(), gs_op_add, this%gs_event)
 
 
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
 
Return the device pointer for an associated Fortran array.
 
Map a Fortran array to a device (allocate and associate)
 
Jacobi preconditioner accelerator backend.
 
subroutine device_jacobi_free(this)
 
subroutine device_jacobi_init(this, coef, dof, gs_h)
 
subroutine device_jacobi_update(this)
 
subroutine device_jacobi_solve(this, z, r, n)
The jacobi preconditioner   where .
 
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
 
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
 
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
 
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
 
Device abstraction, common interface for various accelerators.
 
subroutine, public device_event_sync(event)
Synchronize an event.
 
subroutine, public device_free(x_d)
Deallocate memory on the device.
 
subroutine, public device_event_destroy(event)
Destroy a device event.
 
subroutine, public device_event_create(event, flags)
Create a device event queue.
 
Defines a mapping of the degrees of freedom.
 
integer, parameter, public rp
Global precision used in computations.
 
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
 
Defines a jacobi preconditioner.
 
Defines a canonical Krylov preconditioner.