36  use, 
intrinsic :: iso_c_binding, only : c_int, c_ptr
 
   43          bind(c, name=
'hip_opchsign')
 
   44       use, 
intrinsic :: iso_c_binding
 
   45       type(c_ptr), 
value :: a1_d, a2_d, a3_d
 
   46       integer(c_int) :: gdim, n
 
 
   52          bind(c, name=
'hip_opcolv')
 
   53       use, 
intrinsic :: iso_c_binding
 
   54       type(c_ptr), 
value :: a1_d, a2_d, a3_d, c_d
 
   55       integer(c_int) :: gdim, n
 
 
   60     subroutine hip_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
 
   61          bind(c, name=
'hip_opcolv3c')
 
   62       use, 
intrinsic :: iso_c_binding
 
   64       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
   66       integer(c_int) :: gdim, n
 
 
   71     subroutine hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
 
   72          bind(c, name=
'hip_opadd2cm')
 
   73       use, 
intrinsic :: iso_c_binding
 
   75       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
 
   77       integer(c_int) :: gdim, n
 
 
   82     subroutine hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
 
   83          bind(c, name=
'hip_opadd2col')
 
   84       use, 
intrinsic :: iso_c_binding
 
   85       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
   86       integer(c_int) :: gdim, n
 
 
   92          bind(c, name=
'cuda_opchsign')
 
   93       use, 
intrinsic :: iso_c_binding
 
   94       type(c_ptr), 
value :: a1_d, a2_d, a3_d
 
   95       integer(c_int) :: gdim, n
 
  100     subroutine cuda_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
 
  101          bind(c, name=
'cuda_opcolv')
 
  102       use, 
intrinsic :: iso_c_binding
 
  103       type(c_ptr), 
value :: a1_d, a2_d, a3_d, c_d
 
  104       integer(c_int) :: gdim, n
 
  109     subroutine cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
 
  110          bind(c, name=
'cuda_opcolv3c')
 
  111       use, 
intrinsic :: iso_c_binding
 
  113       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  115       integer(c_int) :: gdim, n
 
  120     subroutine cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
 
  121          bind(c, name=
'cuda_opadd2cm')
 
  122       use, 
intrinsic :: iso_c_binding
 
  124       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
 
  126       integer(c_int) :: gdim, n
 
  131     subroutine cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
 
  132          bind(c, name=
'cuda_opadd2col')
 
  133       use, 
intrinsic :: iso_c_binding
 
  134       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  135       integer(c_int) :: gdim, n
 
  141          bind(c, name=
'opencl_opchsign')
 
  142       use, 
intrinsic :: iso_c_binding
 
  143       type(c_ptr), 
value :: a1_d, a2_d, a3_d
 
  144       integer(c_int) :: gdim, n
 
  150          bind(c, name=
'opencl_opcolv')
 
  151       use, 
intrinsic :: iso_c_binding
 
  152       type(c_ptr), 
value :: a1_d, a2_d, a3_d, c_d
 
  153       integer(c_int) :: gdim, n
 
  158     subroutine opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
 
  159          bind(c, name=
'opencl_opcolv3c')
 
  160       use, 
intrinsic :: iso_c_binding
 
  162       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  164       integer(c_int) :: gdim, n
 
  169     subroutine opencl_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
 
  170          bind(c, name=
'opencl_opadd2cm')
 
  171       use, 
intrinsic :: iso_c_binding
 
  173       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
 
  175       integer(c_int) :: gdim, n
 
  180     subroutine opencl_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
 
  181          bind(c, name=
'opencl_opadd2col')
 
  182       use, 
intrinsic :: iso_c_binding
 
  183       type(c_ptr), 
value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  184       integer(c_int) :: gdim, n
 
  196    type(c_ptr) :: a1_d, a2_d, a3_d
 
  205    call neko_error(
'No device backend configured')
 
 
  211    type(c_ptr) :: a1_d, a2_d, a3_d, c_d
 
  214    call hip_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
 
  220    call neko_error(
'No device backend configured')
 
 
  226                             b1_d, b2_d, b3_d, c_d, d, n, gdim)
 
  227    type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  231    call hip_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
 
  233    call cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
 
  235    call opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
 
  237    call neko_error(
'No device backend configured')
 
 
  243    type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
 
  247    call hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
 
  249    call cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
 
  253    call neko_error(
'No device backend configured')
 
 
  259    type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
 
  262    call hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
 
  264    call cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
 
  268    call neko_error(
'No device backend configured')
 
 
void opencl_opcolv3c(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, real *d, int *gdim, int *n)
 
void opencl_opadd2cm(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, real *c, int *gdim, int *n)
 
void opencl_opcolv(void *a1, void *a2, void *a3, void *c, int *gdim, int *n)
 
void opencl_opchsign(void *a1, void *a2, void *a3, int *gdim, int *n)
 
void opencl_opadd2col(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, int *gdim, int *n)
 
void cuda_opchsign(void *a1, void *a2, void *a3, int *gdim, int *n)
 
void cuda_opcolv3c(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, real *d, int *gdim, int *n)
 
void cuda_opadd2cm(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, real *c, int *gdim, int *n)
 
void cuda_opcolv(void *a1, void *a2, void *a3, void *c, int *gdim, int *n)
 
void cuda_opadd2col(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, int *gdim, int *n)
 
subroutine, public device_opchsign(a1_d, a2_d, a3_d, gdim, n)
 
subroutine, public device_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, n, gdim)
 
subroutine, public device_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
 
subroutine, public device_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, n, gdim)
 
subroutine, public device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
 
integer, parameter, public c_rp
 
integer, parameter, public rp
Global precision used in computations.