34 use,
intrinsic :: iso_c_binding, only: c_ptr, c_int
38 use mpi_f08,
only : mpi_sum, mpi_min, mpi_max, mpi_in_place, mpi_allreduce
77 type(c_ptr) :: a_d, b_d
79 type(c_ptr),
optional :: strm
84 if (
present(strm))
then
99 call neko_error(
'no device backend configured')
106 type(c_ptr) :: a_d, b_d, mask_d
108 type(c_ptr),
optional :: strm
111 if (n .lt. 1 .or. n_mask .lt. 1)
return
113 if (
present(strm))
then
128 call neko_error(
'no device backend configured')
135 type(c_ptr) :: a_d, b_d, mask_d
137 type(c_ptr),
optional :: strm
140 if (n .lt. 1 .or. n_mask .lt. 1)
return
142 if (
present(strm))
then
155 call neko_error(
'no device backend configured')
161 type(c_ptr) :: a_d, b_d, mask_d
163 type(c_ptr),
optional :: strm
166 if (n .lt. 1 .or. n_mask .lt. 1)
return
168 if (
present(strm))
then
183 call neko_error(
'no device backend configured')
189 n2, lx, ly, lz, n_mask, strm)
190 type(c_ptr) :: a_d, b_d, mask_d, facet_d
191 integer :: n1, n2, lx, ly, lz, n_mask
192 type(c_ptr),
optional :: strm
195 if (n_mask .lt. 1)
return
197 if (
present(strm))
then
205 ly, lz, n_mask, strm_)
208 ly, lz, n_mask, strm_)
211 lx, ly, lz, n_mask, strm_)
214 lx, ly, lz, n_mask, strm_)
216 call neko_error(
'no device backend configured')
223 type(c_ptr) :: a_d, b_d, mask_d
225 type(c_ptr),
optional :: strm
228 if (n .lt. 1 .or. n_mask .lt. 1)
return
230 if (
present(strm))
then
245 call neko_error(
'no device backend configured')
251 type(c_ptr) :: a_d, b_d, mask_d
253 type(c_ptr),
optional :: strm
256 if (n .lt. 1 .or. n_mask .lt. 1)
return
258 if (
present(strm))
then
273 call neko_error(
'no device backend configured')
280 type(c_ptr) :: a_d, b_d, mask_d
282 type(c_ptr),
optional :: strm
285 if (n .lt. 1 .or. n_mask .lt. 1)
return
287 if (
present(strm))
then
302 call neko_error(
'no device backend configured')
307 type(c_ptr) :: a_d, b_d, mask_d
309 type(c_ptr),
optional :: strm
312 if (n .lt. 1 .or. n_mask .lt. 1)
return
314 if (
present(strm))
then
325 call neko_error(
'No OpenCL bcknd, masked atomic reduction')
329 call neko_error(
'no device backend configured')
337 real(kind=
rp),
intent(in) :: c
339 type(c_ptr) :: mask_d
341 type(c_ptr),
optional :: strm
344 if (n .lt. 1 .or. n_mask .lt. 1)
return
346 if (
present(strm))
then
361 call neko_error(
'No device backend configured')
369 type(c_ptr),
optional :: strm
374 if (
present(strm))
then
389 call neko_error(
'No device backend configured')
397 type(c_ptr),
optional :: strm
399 real(kind=
rp),
parameter :: one = 1.0_rp
403 if (
present(strm))
then
409#if HAVE_HIP || HAVE_CUDA || HAVE_OPENCL || HAVE_METAL
412 call neko_error(
'No device backend configured')
419 real(kind=
rp),
intent(in) :: c
421 type(c_ptr),
optional :: strm
426 if (
present(strm))
then
441 call neko_error(
'No device backend configured')
447 type(c_ptr) :: a_d, b_d
448 real(kind=
rp),
intent(in) :: c
450 type(c_ptr),
optional :: strm
455 if (
present(strm))
then
470 call neko_error(
'No device backend configured')
477 real(kind=
rp),
intent(in) :: c
479 type(c_ptr),
optional :: strm
482 if (
present(strm))
then
497 call neko_error(
'No device backend configured')
503 type(c_ptr) :: a_d, b_d
504 real(kind=
rp),
intent(in) :: c
506 type(c_ptr),
optional :: strm
509 if (
present(strm))
then
524 call neko_error(
'No device backend configured')
531 real(kind=
rp),
intent(in) :: c
533 type(c_ptr),
optional :: strm
538 if (
present(strm))
then
553 call neko_error(
'No device backend configured')
561 real(kind=
rp),
intent(in) :: c
563 type(c_ptr),
optional :: strm
568 if (
present(strm))
then
583 call neko_error(
'No device backend configured')
590 real(kind=
rp),
intent(in) :: min_val, max_val
592 type(c_ptr),
optional :: strm
595 if (n .lt. 1 .or. max_val .le. min_val)
return
597 if (
present(strm))
then
604 call hip_cwrap(a_d, min_val, max_val, n, strm_)
606 call cuda_cwrap(a_d, min_val, max_val, n, strm_)
612 call neko_error(
'No device backend configured')
619 real(kind=
rp),
intent(in) :: c
621 type(c_ptr),
optional ::strm
626 if (
present(strm))
then
641 call neko_error(
'No device backend configured')
647 type(c_ptr) :: a_d, b_d
649 type(c_ptr),
optional :: strm
654 if (
present(strm))
then
669 call neko_error(
'No device backend configured')
674 type(c_ptr) :: a_d, b_d, c_d, d_d
676 type(c_ptr),
optional :: strm
681 if (
present(strm))
then
688 call hip_add4(a_d, b_d, c_d, d_d, n, strm_)
690 call cuda_add4(a_d, b_d, c_d, d_d, n, strm_)
696 call neko_error(
'No device backend configured')
701 type(c_ptr) :: a_d, b_d
704 type(c_ptr),
optional :: strm
709 if (
present(strm))
then
724 call neko_error(
'No device backend configured')
731 type(c_ptr) :: a_d, b_d
734 type(c_ptr),
optional :: strm
739 if (
present(strm))
then
754 call neko_error(
'No device backend configured')
760 type(c_ptr) :: a_d, b_d
763 type(c_ptr),
optional :: strm
768 if (
present(strm))
then
783 call neko_error(
'No device backend configured')
789 type(c_ptr) :: a_d, b_d, c_d
791 type(c_ptr),
optional :: strm
796 if (
present(strm))
then
803 call hip_add3(a_d, b_d, c_d, n, strm_)
811 call neko_error(
'No device backend configured')
817 type(c_ptr) :: a_d, b_d, c_d
818 real(kind=
rp) :: c1, c2
820 type(c_ptr),
optional :: strm
825 if (
present(strm))
then
832 call hip_add3s2(a_d, b_d, c_d, c1, c2, n, strm_)
840 call neko_error(
'No device backend configured')
846 type(c_ptr) :: a_d, b_d, c_d, d_d
847 real(kind=
rp) :: c1, c2, c3
849 type(c_ptr),
optional :: strm
854 if (
present(strm))
then
861 call hip_add4s3(a_d, b_d, c_d, d_d, c1, c2, c3, n, strm_)
863 call cuda_add4s3(a_d, b_d, c_d, d_d, c1, c2, c3, n, strm_)
867 call metal_add4s3(a_d, b_d, c_d, d_d, c1, c2, c3, n, strm_)
869 call neko_error(
'No device backend configured')
874 subroutine device_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2 , c3, c4, n, strm)
875 type(c_ptr) :: a_d, b_d, c_d, d_d, e_d
876 real(kind=
rp) :: c1, c2, c3, c4
878 type(c_ptr),
optional :: strm
883 if (
present(strm))
then
890 call hip_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm_)
892 call cuda_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm_)
894 call opencl_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm_)
896 call metal_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm_)
898 call neko_error(
'No device backend configured')
906 type(c_ptr),
optional :: strm
911 if (
present(strm))
then
926 call neko_error(
'No device backend configured')
932 type(c_ptr) :: a_d, b_d
934 type(c_ptr),
optional :: strm
939 if (
present(strm))
then
954 call neko_error(
'No device backend configured')
960 type(c_ptr) :: a_d, b_d, c_d
962 type(c_ptr),
optional :: strm
965 if (
present(strm))
then
977 call neko_error(
'opencl_invcol3 not implemented')
981 call neko_error(
'No device backend configured')
987 type(c_ptr) :: a_d, b_d
989 type(c_ptr),
optional :: strm
992 if (
present(strm))
then
1008 call neko_error(
'No device backend configured')
1014 type(c_ptr) :: a_d, b_d, c_d
1016 type(c_ptr),
optional :: strm
1017 type(c_ptr) :: strm_
1019 if (n .lt. 1)
return
1021 if (
present(strm))
then
1028 call hip_col3(a_d, b_d, c_d, n, strm_)
1036 call neko_error(
'No device backend configured')
1042 type(c_ptr) :: a_d, b_d, c_d
1044 type(c_ptr),
optional :: strm
1045 type(c_ptr) :: strm_
1047 if (n .lt. 1)
return
1049 if (
present(strm))
then
1064 call neko_error(
'No device backend configured')
1070 type(c_ptr) :: a_d, b_d
1072 type(c_ptr),
optional :: strm
1073 type(c_ptr) :: strm_
1075 if (n .lt. 1)
return
1077 if (
present(strm))
then
1092 call neko_error(
'No device backend configured')
1098 type(c_ptr) :: a_d, b_d, c_d
1100 type(c_ptr),
optional :: strm
1101 type(c_ptr) :: strm_
1103 if (n .lt. 1)
return
1105 if (
present(strm))
then
1112 call hip_sub3(a_d, b_d, c_d, n, strm_)
1120 call neko_error(
'No device backend configured')
1126 type(c_ptr) :: a_d, b_d, c_d
1128 type(c_ptr),
optional :: strm
1129 type(c_ptr) :: strm_
1131 if (n .lt. 1)
return
1133 if (
present(strm))
then
1148 call neko_error(
'No device backend configured')
1154 type(c_ptr) :: a_d, b_d, c_d, d_d
1156 type(c_ptr),
optional :: strm
1157 type(c_ptr) :: strm_
1159 if (n .lt. 1)
return
1161 if (
present(strm))
then
1176 call neko_error(
'No device backend configured')
1182 type(c_ptr) :: a_d, b_d, c_d
1185 type(c_ptr),
optional :: strm
1186 type(c_ptr) :: strm_
1188 if (n .lt. 1)
return
1190 if (
present(strm))
then
1205 call neko_error(
'No device backend configured')
1211 subroutine device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
1212 type(c_ptr) :: dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d
1214 type(c_ptr),
optional :: strm
1215 type(c_ptr) :: strm_
1217 if (n .lt. 1)
return
1219 if (
present(strm))
then
1226 call hip_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm_)
1228 call cuda_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm_)
1230 call opencl_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm_)
1232 call metal_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm_)
1234 call neko_error(
'No device backend configured')
1241 w1_d, w2_d, w3_d, n, strm)
1242 type(c_ptr) :: u1_d, u2_d, u3_d
1243 type(c_ptr) :: v1_d, v2_d, v3_d
1244 type(c_ptr) :: w1_d, w2_d, w3_d
1246 type(c_ptr),
optional :: strm
1247 type(c_ptr) :: strm_
1249 if (n .lt. 1)
return
1251 if (
present(strm))
then
1258 call hip_vcross(u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, &
1259 w1_d, w2_d, w3_d, n, strm_)
1261 call cuda_vcross(u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, &
1262 w1_d, w2_d, w3_d, n, strm_)
1265 w1_d, w2_d, w3_d, n, strm_)
1267 call metal_vcross(u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, &
1268 w1_d, w2_d, w3_d, n, strm_)
1270 call neko_error(
'No device backend configured')
1277 type(c_ptr) :: u_d, v_d, w_d
1279 type(c_ptr),
optional :: strm
1280 type(c_ptr) :: strm_
1281 real(kind=
rp) :: res
1285 if (n .lt. 1)
return
1287 if (
present(strm))
then
1294 res =
hip_vlsc3(u_d, v_d, w_d, n, strm_)
1304 call neko_error(
'No device backend configured')
1310 type(c_ptr) :: a_d, b_d, c_d
1312 type(c_ptr),
optional :: strm
1313 type(c_ptr) :: strm_
1314 real(kind=
rp) :: res
1316 if (
present(strm))
then
1324 res =
hip_glsc3(a_d, b_d, c_d, n, strm_)
1332 call neko_error(
'No device backend configured')
1335#ifndef HAVE_DEVICE_MPI
1337 call mpi_allreduce(mpi_in_place, res, 1, &
1344 type(c_ptr),
value :: w_d, v_d_d, mult_d
1345 integer(c_int) :: j, n
1347 type(c_ptr),
optional :: strm
1348 type(c_ptr) :: strm_
1351 if (
present(strm))
then
1366 call neko_error(
'No device backend configured')
1369#ifndef HAVE_DEVICE_MPI
1371 call mpi_allreduce(mpi_in_place, h, j, &
1378 type(c_ptr),
value :: y_d, x_d_d, a_d
1379 integer(c_int) :: j, n
1380 type(c_ptr),
optional :: strm
1381 type(c_ptr) :: strm_
1383 if (n .lt. 1)
return
1385 if (
present(strm))
then
1400 call neko_error(
'No device backend configured')
1406 type(c_ptr) :: a_d, b_d
1408 real(kind=
rp) :: res
1409 type(c_ptr),
optional :: strm
1410 type(c_ptr) :: strm_
1412 if (
present(strm))
then
1428 call neko_error(
'No device backend configured')
1431#ifndef HAVE_DEVICE_MPI
1433 call mpi_allreduce(mpi_in_place, res, 1, &
1442 type(c_ptr),
intent(in) :: a_d, b_d
1443 integer,
intent(in) :: n
1445 real(kind=
rp) :: res
1446 type(c_ptr),
optional :: strm
1447 type(c_ptr) :: strm_
1449 if (
present(strm))
then
1465 call neko_error(
'No device backend configured')
1468#ifndef HAVE_DEVICE_MPI
1470 call mpi_allreduce(mpi_in_place, res, 1, &
1482 real(kind=
rp) :: res
1483 type(c_ptr),
optional :: strm
1484 type(c_ptr) :: strm_
1486 if (
present(strm))
then
1502 call neko_error(
'No device backend configured')
1505#ifndef HAVE_DEVICE_MPI
1507 call mpi_allreduce(mpi_in_place, res, 1, &
1517 real(kind=
rp) :: res, ninf
1518 type(c_ptr),
optional :: strm
1519 type(c_ptr) :: strm_
1526 if (
present(strm))
then
1532 ninf = -huge(0.0_rp)
1542 call neko_error(
'No device backend configured')
1545#ifndef HAVE_DEVICE_MPI
1547 call mpi_allreduce(mpi_in_place, res, 1, &
1557 real(kind=
rp) :: res, pinf
1558 type(c_ptr),
optional :: strm
1559 type(c_ptr) :: strm_
1566 if (
present(strm))
then
1582 call neko_error(
'No device backend configured')
1585#ifndef HAVE_DEVICE_MPI
1587 call mpi_allreduce(mpi_in_place, res, 1, &
1594 integer,
intent(in) :: n
1596 type(c_ptr),
optional :: strm
1597 type(c_ptr) :: strm_
1599 if (n .lt. 1)
return
1601 if (
present(strm))
then
1616 call neko_error(
'No device backend configured')
1627 type(c_ptr) :: a_d, b_d
1629 type(c_ptr),
optional :: strm
1630 type(c_ptr) :: strm_
1632 if (n .lt. 1)
return
1634 if (
present(strm))
then
1649 call neko_error(
'No device backend configured')
1656 type(c_ptr) :: a_d, b_d, c_d
1658 type(c_ptr),
optional :: strm
1659 type(c_ptr) :: strm_
1661 if (n .lt. 1)
return
1663 if (
present(strm))
then
1678 call neko_error(
'No device backend configured')
1687 real(kind=
rp),
intent(in) :: c
1689 type(c_ptr),
optional :: strm
1690 type(c_ptr) :: strm_
1692 if (n .lt. 1)
return
1694 if (
present(strm))
then
1709 call neko_error(
'No device backend configured')
1717 type(c_ptr) :: a_d, b_d
1718 real(kind=
rp),
intent(in) :: c
1720 type(c_ptr),
optional :: strm
1721 type(c_ptr) :: strm_
1723 if (n .lt. 1)
return
1725 if (
present(strm))
then
1740 call neko_error(
'No device backend configured')
1751 type(c_ptr) :: a_d, b_d
1753 type(c_ptr),
optional :: strm
1754 type(c_ptr) :: strm_
1756 if (n .lt. 1)
return
1758 if (
present(strm))
then
1773 call neko_error(
'No device backend configured')
1780 type(c_ptr) :: a_d, b_d, c_d
1782 type(c_ptr),
optional :: strm
1783 type(c_ptr) :: strm_
1785 if (n .lt. 1)
return
1787 if (
present(strm))
then
1802 call neko_error(
'No device backend configured')
1811 real(kind=
rp),
intent(in) :: c
1813 type(c_ptr),
optional :: strm
1814 type(c_ptr) :: strm_
1816 if (n .lt. 1)
return
1818 if (
present(strm))
then
1833 call neko_error(
'No device backend configured')
1841 type(c_ptr) :: a_d, b_d
1842 real(kind=
rp),
intent(in) :: c
1844 type(c_ptr),
optional :: strm
1845 type(c_ptr) :: strm_
1847 if (n .lt. 1)
return
1849 if (
present(strm))
then
1864 call neko_error(
'No device backend configured')
1874 type(c_ptr),
intent(inout) :: a_d
1875 integer,
intent(in) :: c
1876 integer,
intent(in) :: n
1877 type(c_ptr),
optional :: strm
1878 type(c_ptr) :: strm_
1879 if (n .lt. 1)
return
1881 if (
present(strm))
then
1896 call neko_error(
'No device backend configured')
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
integer, public pe_size
MPI size of communicator.
type(mpi_comm), public neko_comm
MPI communicator.
subroutine, public device_pwmin2(a_d, b_d, n, strm)
Compute the point-wise minimum of two vectors .
subroutine, public device_masked_scatter_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_add2s1(a_d, b_d, c1, n, strm)
subroutine, public device_add2s2_many(y_d, x_d_d, a_d, j, n, strm)
subroutine, public device_add4s3(a_d, b_d, c_d, d_d, c1, c2, c3, n, strm)
Returns .
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_masked_scatter_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
real(kind=rp) function, public device_glmax(a_d, n, strm)
Max of a vector of length n.
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_pwmax3(a_d, b_d, c_d, n, strm)
Compute the point-wise maximum of two vectors .
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n, strm)
Returns .
subroutine, public device_masked_atomic_reduction_0(a_d, b_d, mask_d, n, n_mask, strm)
subroutine, public device_cpwmax3(a_d, b_d, c, n, strm)
Compute the point-wise maximum of a vector and a scalar .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_rone(a_d, n, strm)
Set all elements to one.
subroutine, public device_vcross(u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, w1_d, w2_d, w3_d, n, strm)
Compute a cross product (3-d version) assuming vector components etc.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
real(kind=rp) function, public device_glsubnorm(a_d, b_d, n, strm)
Returns the norm of the difference of two vectors .
subroutine, public device_glsc3_many(h, w_d, v_d_d, mult_d, j, n, strm)
subroutine, public device_cfill_mask(a_d, c, n, mask_d, n_mask, strm)
Fill a constant to a masked vector. .
subroutine, public device_cadd2(a_d, b_d, c, n, strm)
Add a scalar to vector .
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_invcol3(a_d, b_d, c_d, n, strm)
Vector division .
subroutine, public device_pwmin3(a_d, b_d, c_d, n, strm)
Compute the point-wise minimum of two vectors .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_masked_gather_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
real(kind=rp) function, public device_vlsc3(u_d, v_d, w_d, n, strm)
Compute multiplication sum .
subroutine, public device_cdiv2(a_d, b_d, c, n, strm)
Division of constant c by array .
subroutine, public device_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm)
Returns .
subroutine, public device_add4(a_d, b_d, c_d, d_d, n, strm)
subroutine, public device_cpwmax2(a_d, c, n, strm)
Compute the point-wise maximum of a vector and a scalar .
subroutine, public device_masked_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Copy a masked vector .
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_cdiv(a_d, c, n, strm)
Division of constant c by array .
subroutine, public device_absval(a_d, n, strm)
subroutine device_iadd(a_d, c, n, strm)
Add an integer scalar to vector .
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n, strm)
Returns .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n, strm)
Weighted inner product .
subroutine, public device_cpwmin3(a_d, b_d, c, n, strm)
Compute the point-wise minimum of a vector and a scalar .
subroutine, public device_cwrap(a_d, min_val, max_val, n, strm)
Wrap value around a range (min, max)
subroutine, public device_face_masked_gather_copy_0(a_d, b_d, mask_d, facet_d, n1, n2, lx, ly, lz, n_mask, strm)
Gather a face-local SEM field .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_masked_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Copy a masked vector .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
subroutine, public device_addcol4(a_d, b_d, c_d, d_d, n, strm)
Returns .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_cpwmin2(a_d, c, n, strm)
Compute the point-wise minimum of a vector and a scalar .
subroutine, public device_add3(a_d, b_d, c_d, n, strm)
Vector addition .
real(kind=rp) function, public device_glmin(a_d, n, strm)
Min of a vector of length n.
subroutine, public device_addcol3s2(a_d, b_d, c_d, s, n, strm)
Returns .
subroutine, public device_pwmax2(a_d, b_d, n, strm)
Compute the point-wise maximum of two vectors .
subroutine device_radd(a_d, c, n, strm)
Add a scalar to vector .
Device abstraction, common interface for various accelerators.
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.