63 use mpi_f08,
only: mpi_min, mpi_max, mpi_sum, mpi_in_place, mpi_integer, &
69 real(kind=
rp),
public,
parameter ::
neko_eps = epsilon(1.0_rp)
75 real(kind=
rp),
public,
parameter ::
pi = 4._rp*atan(1._rp)
111 vcross,
vdot2,
vdot3,
vlsc3,
vlsc2,
add2,
add3,
add4,
sub2,
sub3, &
123 real(kind=
sp),
intent(in) :: x
124 real(kind=
sp),
intent(in) :: y
125 real(kind=
sp),
intent(in),
optional :: tol
128 if (
present(tol))
then
138 real(kind=
dp),
intent(in) :: x
139 real(kind=
dp),
intent(in) :: y
140 real(kind=
dp),
intent(in),
optional :: tol
143 if (
present(tol))
then
153 real(kind=
qp),
intent(in) :: x
154 real(kind=
qp),
intent(in) :: y
155 real(kind=
qp),
intent(in),
optional :: tol
158 if (
present(tol))
then
168 real(kind=
sp),
intent(in) :: x
169 real(kind=
sp),
intent(in) :: y
170 real(kind=
sp),
intent(in),
optional :: eps
172 if (
present(eps))
then
173 srelcmp = abs(x - y) .le. eps*abs(y)
182 real(kind=
dp),
intent(in) :: x
183 real(kind=
dp),
intent(in) :: y
184 real(kind=
dp),
intent(in),
optional :: eps
186 if (
present(eps))
then
187 drelcmp = abs(x - y) .le. eps*abs(y)
197 real(kind=
qp),
intent(in) :: x
198 real(kind=
qp),
intent(in) :: y
199 real(kind=
qp),
intent(in),
optional :: eps
201 if (
present(eps))
then
202 qrelcmp = abs(x - y)/abs(y) .lt. eps
211 integer,
intent(in) :: n
212 real(kind=
rp),
dimension(n),
intent(inout) :: a
222 integer,
intent(in) :: n
223 integer,
dimension(n),
intent(inout) :: a
233 integer,
intent(in) :: m, n, e
234 real(kind=
rp),
intent(inout) :: a(m,n)
244 integer,
intent(in) :: n
245 real(kind=
rp),
dimension(n),
intent(inout) :: a
255 integer,
intent(in) :: n
256 real(kind=
rp),
dimension(n),
intent(in) :: b
257 real(kind=
rp),
dimension(n),
intent(inout) :: a
274 integer,
intent(in) :: n, n_mask
275 real(kind=
rp),
dimension(n),
intent(in) :: b
276 real(kind=
rp),
dimension(n),
intent(inout) :: a
277 integer,
dimension(0:n_mask) ::
mask
295 integer,
intent(in) :: n, n_mask
296 real(kind=
rp),
dimension(n),
intent(in) :: b
297 real(kind=
rp),
dimension(n),
intent(inout) :: a
298 integer,
dimension(n_mask) ::
mask
318 integer,
intent(in) :: n, n_mask
319 real(kind=
rp),
dimension(n),
intent(in) :: b
320 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
321 integer,
dimension(0:n_mask) ::
mask
341 integer,
intent(in) :: n, n_mask
342 real(kind=
rp),
dimension(n),
intent(in) :: b
343 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
344 integer,
dimension(n_mask) ::
mask
364 integer,
intent(in) :: n, n_mask
365 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
366 real(kind=
rp),
dimension(n),
intent(inout) :: a
367 integer,
dimension(0:n_mask) ::
mask
387 integer,
intent(in) :: n, n_mask
388 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
389 real(kind=
rp),
dimension(n),
intent(inout) :: a
390 integer,
dimension(n_mask) ::
mask
403 integer,
intent(in) :: n, n_mask
404 real(kind=
rp),
dimension(n),
intent(inout) :: a
405 real(kind=
rp),
intent(in) :: c
406 integer,
dimension(n_mask),
intent(in) ::
mask
417 integer,
intent(in) :: n
418 real(kind=
rp),
dimension(n),
intent(inout) :: a
419 real(kind=
rp),
intent(in) :: c
429 integer,
intent(in) :: n
430 real(kind=
rp),
dimension(n),
intent(inout) :: a
431 real(kind=
rp),
dimension(n),
intent(in) :: b
432 real(kind=
rp),
intent(in) :: c
443 integer,
intent(in) :: n
444 real(kind=
rp),
dimension(n),
intent(inout) :: a
445 real(kind=
rp),
intent(in) :: c
455 integer,
intent(in) :: n
456 real(kind=
rp),
dimension(n),
intent(inout) :: a
457 real(kind=
rp),
dimension(n),
intent(in) :: b
458 real(kind=
rp),
intent(in) :: c
468 integer,
intent(in) :: n
469 real(kind=
rp),
dimension(n),
intent(inout) :: a
470 real(kind=
rp),
intent(in) :: s
480 integer,
intent(in) :: n
481 real(kind=
rp),
dimension(n),
intent(inout) :: a
482 real(kind=
rp),
dimension(n),
intent(in) :: b
483 real(kind=
rp),
intent(in) :: s
493 integer,
intent(in) :: n
494 real(kind=
rp),
dimension(n),
intent(inout) :: a
495 real(kind=
rp),
intent(in) :: c
505 integer,
intent(in) :: n
506 real(kind=
rp),
dimension(n) :: a
514 call mpi_allreduce(mpi_in_place, tmp, 1, &
522 integer,
intent(in) :: n
523 real(kind=
rp),
dimension(n) :: a
531 call mpi_allreduce(tmp,
glmax, 1, &
537 integer,
intent(in) :: n
538 integer,
dimension(n) :: a
546 call mpi_allreduce(tmp,
glimax, 1, &
552 integer,
intent(in) :: n
553 real(kind=
rp),
dimension(n) :: a
561 call mpi_allreduce(tmp,
glmin, 1, &
567 integer,
intent(in) :: n
568 integer,
dimension(n) :: a
576 call mpi_allreduce(tmp,
glimin, 1, &
585 integer,
intent(in) :: n
586 real(kind=
rp),
dimension(n),
intent(inout) :: a
598 real(kind=
rp),
intent(in) :: vec(n)
599 real(kind=
rp) :: tmax
602 tmax =
max(tmax, vec(i))
608 integer,
intent(in) :: n
609 real(kind=
rp),
intent(in) :: vec(n)
610 real(kind=
rp) :: tmin
614 tmin = min(tmin, vec(i))
620 integer,
intent(in) :: n
621 real(kind=
rp),
dimension(n),
intent(inout) :: a
625 a(i) = 1.0_xp /
real(a(i),
xp)
632 integer,
intent(in) :: n
633 real(kind=
rp),
dimension(n),
intent(inout) :: a
634 real(kind=
rp),
dimension(n),
intent(in) :: b,c
638 a(i) =
real(b(i),
xp) / c(i)
645 integer,
intent(in) :: n
646 real(kind=
rp),
dimension(n),
intent(inout) :: a
647 real(kind=
rp),
dimension(n),
intent(in) :: b
651 a(i) = 1.0_xp /
real(b(i),
xp)
658 subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
659 integer,
intent(in) :: n
660 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
661 real(kind=
rp),
dimension(n),
intent(in) :: w1, w2, w3
662 real(kind=
rp),
dimension(n),
intent(out) :: u1, u2, u3
666 u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
667 u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
668 u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
675 subroutine vdot2(dot, u1, u2, v1, v2, n)
676 integer,
intent(in) :: n
677 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2
678 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2
679 real(kind=
rp),
dimension(n),
intent(out) :: dot
682 dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
689 subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
690 integer,
intent(in) :: n
691 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2, u3
692 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
693 real(kind=
rp),
dimension(n),
intent(out) :: dot
697 dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
703 function vlsc3(u, v, w, n)
result(s)
704 integer,
intent(in) :: n
705 real(kind=
rp),
dimension(n),
intent(in) :: u, v, w
711 s = s + u(i)*v(i)*w(i)
718 integer,
intent(in) :: n
719 real(kind=
rp),
dimension(n),
intent(in) :: u, v
732 integer,
intent(in) :: n
733 real(kind=
rp),
dimension(n),
intent(inout) :: a
734 real(kind=
rp),
dimension(n),
intent(in) :: b
745 integer,
intent(in) :: n
746 real(kind=
rp),
dimension(n),
intent(inout) :: a
747 real(kind=
rp),
dimension(n),
intent(in) :: b
748 real(kind=
rp),
dimension(n),
intent(in) :: c
759 integer,
intent(in) :: n
760 real(kind=
rp),
dimension(n),
intent(out) :: a
761 real(kind=
rp),
dimension(n),
intent(in) :: d
762 real(kind=
rp),
dimension(n),
intent(in) :: c
763 real(kind=
rp),
dimension(n),
intent(in) :: b
767 a(i) = b(i) + c(i) + d(i)
774 integer,
intent(in) :: n
775 real(kind=
rp),
dimension(n),
intent(inout) :: a
776 real(kind=
rp),
dimension(n),
intent(inout) :: b
787 integer,
intent(in) :: n
788 real(kind=
rp),
dimension(n),
intent(inout) :: a
789 real(kind=
rp),
dimension(n),
intent(in) :: b
790 real(kind=
rp),
dimension(n),
intent(in) :: c
803 integer,
intent(in) :: n
804 real(kind=
rp),
dimension(n),
intent(inout) :: a
805 real(kind=
rp),
dimension(n),
intent(inout) :: b
806 real(kind=
rp),
intent(in) :: c1
810 a(i) = c1 * a(i) + b(i)
818 integer,
intent(in) :: n
819 real(kind=
rp),
dimension(n),
intent(inout) :: a
820 real(kind=
rp),
dimension(n),
intent(in) :: b
821 real(kind=
rp),
intent(in) :: c1
825 a(i) = a(i) + c1 * b(i)
832 integer,
intent(in) :: n
833 real(kind=
rp),
dimension(n),
intent(inout) :: a
834 real(kind=
rp),
dimension(n),
intent(in) :: b
835 real(kind=
rp),
intent(in) :: c1
839 a(i) = a(i) + c1 * ( b(i) * b(i) )
846 integer,
intent(in) :: n
847 real(kind=
rp),
dimension(n),
intent(inout) :: a
848 real(kind=
rp),
dimension(n),
intent(in) :: b
852 a(i) =
real(a(i),
xp) /b(i)
860 integer,
intent(in) :: n
861 real(kind=
rp),
dimension(n),
intent(inout) :: a
862 real(kind=
rp),
dimension(n),
intent(in) :: b
873 integer,
intent(in) :: n
874 real(kind=
rp),
dimension(n),
intent(inout) :: a
875 real(kind=
rp),
dimension(n),
intent(in) :: b
876 real(kind=
rp),
dimension(n),
intent(in) :: c
887 integer,
intent(in) :: n
888 real(kind=
rp),
dimension(n),
intent(inout) :: a
889 real(kind=
rp),
dimension(n),
intent(in) :: b
890 real(kind=
rp),
dimension(n),
intent(in) :: c
894 a(i) = a(i) - b(i) * c(i)
901 integer,
intent(in) :: n
902 real(kind=
rp),
dimension(n),
intent(inout) :: a
903 real(kind=
rp),
dimension(n),
intent(in) :: b
904 real(kind=
rp),
dimension(n),
intent(in) :: c
905 real(kind=
rp),
intent(in) :: c1, c2
909 a(i) = c1 * b(i) + c2 * c(i)
917 integer,
intent(in) :: n
918 real(kind=
rp),
dimension(n),
intent(inout) :: a
919 real(kind=
rp),
dimension(n),
intent(in) :: b
920 real(kind=
rp),
dimension(n),
intent(in) :: c
921 real(kind=
rp),
dimension(n),
intent(in) :: d
925 a(i) = a(i) - b(i) * c(i) * d(i)
932 integer,
intent(in) :: n
933 real(kind=
rp),
dimension(n),
intent(inout) :: a
934 real(kind=
rp),
dimension(n),
intent(in) :: b
935 real(kind=
rp),
dimension(n),
intent(in) :: c
939 a(i) = a(i) + b(i) * c(i)
946 integer,
intent(in) :: n
947 real(kind=
rp),
dimension(n),
intent(inout) :: a
948 real(kind=
rp),
dimension(n),
intent(in) :: b
949 real(kind=
rp),
dimension(n),
intent(in) :: c
950 real(kind=
rp),
dimension(n),
intent(in) :: d
954 a(i) = a(i) + b(i) * c(i) * d(i)
961 integer,
intent(in) :: n
962 real(kind=
rp),
dimension(n),
intent(inout) :: a
963 real(kind=
rp),
dimension(n),
intent(in) :: b
964 real(kind=
rp),
dimension(n),
intent(in) :: c
965 real(kind=
rp),
dimension(n),
intent(in) :: d
966 real(kind=
rp),
dimension(n),
intent(in) :: e
970 a(i) = b(i)*c(i)-d(i)*e(i)
977 integer,
intent(in) :: n
978 real(kind=
rp),
dimension(n),
intent(inout) :: a
979 real(kind=
rp),
dimension(n),
intent(in) :: b
980 real(kind=
rp),
dimension(n),
intent(in) :: c
981 real(kind=
rp),
intent(in) :: c1, c2
985 a(i) = b(i) + c1*(a(i)-c2*c(i))
992 integer,
intent(in) :: n
993 real(kind=
rp),
dimension(n),
intent(inout) :: a
994 real(kind=
rp),
dimension(n),
intent(in) :: b
995 real(kind=
rp),
dimension(n),
intent(in) :: c
996 real(kind=
rp),
intent(in) :: c1, c2
1000 a(i) = a(i) + c1*b(i)+c2*c(i)
1007 integer,
intent(in) :: n
1008 real(kind=
rp),
dimension(n),
intent(in) :: a
1009 real(kind=
rp),
dimension(n),
intent(in) :: b
1011 real(kind=
xp) :: tmp
1016 tmp = tmp + a(i) * b(i)
1019 call mpi_allreduce(mpi_in_place, tmp, 1, &
1026 integer,
intent(in) :: n
1027 real(kind=
rp),
dimension(n),
intent(in) :: a
1028 real(kind=
rp),
dimension(n),
intent(in) :: b
1029 real(kind=
rp),
dimension(n),
intent(in) :: c
1031 real(kind=
xp) :: tmp
1036 tmp = tmp + a(i) * b(i) * c(i)
1039 call mpi_allreduce(mpi_in_place, tmp, 1, &
1045 integer,
intent(in) :: n
1046 real(kind=
rp),
dimension(n),
intent(in) :: a
1047 real(kind=
rp),
dimension(n),
intent(in) :: b
1048 real(kind=
rp),
dimension(n),
intent(in) :: c
1049 real(kind=
rp),
dimension(n),
intent(in) :: d
1051 real(kind=
xp) :: tmp
1056 tmp = tmp + a(i) * b(i) * c(i) * d(i)
1059 call mpi_allreduce(mpi_in_place, tmp, 1, &
1068 integer,
intent(in) :: n
1069 real(kind=
rp),
dimension(n),
intent(in) :: a
1070 real(kind=
rp),
dimension(n),
intent(in) :: b
1072 real(kind=
xp) :: tmp
1077 tmp = tmp + (a(i) - b(i))**2
1080 call mpi_allreduce(mpi_in_place, tmp, 1, &
1092 integer,
intent(in) :: n
1093 real(kind=
rp),
intent(inout) :: a(n)
1094 integer,
intent(out) :: ind(n)
1096 integer :: j, ir, i, ii, l
1102 if (n .le. 1)
return
1126 do while (j .le. ir)
1128 if ( a(j) .lt. a(j+1) ) j = j + 1
1130 if (aa .lt. a(j))
then
1150 integer,
intent(in) :: n
1151 integer(i4),
intent(inout) :: a(n)
1152 integer,
intent(out) :: ind(n)
1154 integer :: j, ir, i, ii, l
1160 if (n .le. 1)
return
1183 do while (j .le. ir)
1185 if ( a(j) .lt. a(j + 1) ) j = j + 1
1187 if (aa .lt. a(j))
then
1206 integer,
intent(in) :: n
1207 real(kind=
rp),
intent(inout) :: b(n)
1208 integer,
intent(in) :: ind(n)
1209 real(kind=
rp) :: temp(n)
1226 integer,
intent(in) :: n
1227 integer(i4),
intent(inout) :: b(n)
1228 integer,
intent(in) :: ind(n)
1229 integer(i4) :: temp(n)
1246 integer,
intent(in) :: n
1247 real(kind=
rp),
intent(inout) :: b(n)
1248 integer,
intent(in) :: ind(n)
1249 real(kind=
rp) :: temp(n)
1266 integer,
intent(in) :: n
1267 integer(i4),
intent(inout) :: b(n)
1268 integer,
intent(in) :: ind(n)
1269 integer(i4) :: temp(n)
1286 integer,
intent(in) :: n
1287 real(kind=
rp),
intent(inout) :: b(n)
1288 integer,
intent(inout) :: ind(n)
1289 real(kind=
rp) :: temp(n)
1290 integer :: tempind(n)
1296 tempind(jj) = ind(i)
1309 integer,
intent(in) :: n
1310 integer(i4),
intent(inout) :: b(n)
1311 integer,
intent(inout) :: ind(n)
1312 integer(i4) :: temp(n)
1313 integer :: tempind(n)
1319 tempind(jj) = ind(i)
1331 integer,
intent(in) :: n
1332 real(kind=
rp),
dimension(n),
intent(inout) :: a
1344 integer,
intent(in) :: n
1345 real(kind=
rp),
dimension(n),
intent(inout) :: a
1346 real(kind=
rp),
dimension(n),
intent(in) :: b
1350 a(i) =
max(a(i), b(i))
1356 integer,
intent(in) :: n
1357 real(kind=
rp),
dimension(n),
intent(inout) :: a
1358 real(kind=
rp),
dimension(n),
intent(in) :: b, c
1362 a(i) =
max(b(i), c(i))
1368 integer,
intent(in) :: n
1369 real(kind=
rp),
dimension(n),
intent(inout) :: a
1370 real(kind=
rp),
intent(in) :: b
1380 integer,
intent(in) :: n
1381 real(kind=
rp),
dimension(n),
intent(inout) :: a
1382 real(kind=
rp),
dimension(n),
intent(in) :: b
1383 real(kind=
rp),
intent(in) :: c
1393 integer,
intent(in) :: n
1394 real(kind=
rp),
dimension(n),
intent(inout) :: a
1395 real(kind=
rp),
dimension(n),
intent(in) :: b
1399 a(i) = min(a(i), b(i))
1405 integer,
intent(in) :: n
1406 real(kind=
rp),
dimension(n),
intent(inout) :: a
1407 real(kind=
rp),
dimension(n),
intent(in) :: b, c
1411 a(i) = min(b(i), c(i))
1417 integer,
intent(in) :: n
1418 real(kind=
rp),
dimension(n),
intent(inout) :: a
1419 real(kind=
rp),
intent(in) :: b
1429 integer,
intent(in) :: n
1430 real(kind=
rp),
dimension(n),
intent(inout) :: a
1431 real(kind=
rp),
dimension(n),
intent(in) :: b
1432 real(kind=
rp),
intent(in) :: c
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
type(mpi_comm), public neko_comm
MPI communicator.
type(mpi_datatype), public mpi_extra_precision
Object for handling masks in Neko.
subroutine, public cmult(a, c, n)
Multiplication by constant c .
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
subroutine, public invcol2(a, b, n)
Vector division .
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
pure logical function sabscmp(x, y, tol)
Return single precision absolute comparison .
subroutine pwmin_sca3(a, b, c, n)
Point-wise minimum of scalar and vector .
real(kind=rp), parameter, public pi
pure logical function qabscmp(x, y, tol)
Return double precision absolute comparison .
subroutine pwmax_scal2(a, b, n)
Point-wise maximum of scalar and vector .
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
subroutine, public masked_scatter_copy(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
subroutine, public invers2(a, b, n)
Compute inverted vector .
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
subroutine pwmin_sca2(a, b, n)
Point-wise minimum of scalar and vector .
subroutine, public cadd(a, s, n)
Add a scalar to vector .
subroutine, public masked_copy(a, b, mask, n, n_mask)
Copy a masked vector .
subroutine reorddp(b, ind, n)
reorder double precision array - inverse of swap
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
real(kind=rp) function, public glsc4(a, b, c, d, n)
subroutine, public cdiv2(a, b, c, n)
Division of constant c by elements of a .
subroutine pwmax_vec3(a, b, c, n)
Point-wise maximum of two vectors .
subroutine swapdp(b, ind, n)
sort double precision array acording to ind vector
subroutine flipvi4(b, ind, n)
Flip single integer vector b and ind.
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
subroutine, public masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
subroutine, public subcol3(a, b, c, n)
Returns .
subroutine, public rone(a, n)
Set all elements to one.
subroutine flipvdp(b, ind, n)
Flip double precision vector b and ind.
subroutine, public masked_gather_copy(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
subroutine, public add3(a, b, c, n)
Vector addition .
subroutine swapi4(b, ind, n)
sort single integer array acording to ind vector
integer function, public glimin(a, n)
Min of an integer vector of length n.
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
subroutine, public sub3(a, b, c, n)
Vector subtraction .
subroutine, public addcol4(a, b, c, d, n)
Returns .
subroutine, public add2(a, b, n)
Vector addition .
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
subroutine, public absval(a, n)
Take the absolute value of an array.
subroutine, public invcol3(a, b, c, n)
Invert a vector .
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
subroutine, public subcol4(a, b, c, d, n)
Returns .
subroutine sorti4(a, ind, n)
Heap Sort for single integer arrays.
subroutine, public addcol3(a, b, c, n)
Returns .
subroutine, public invcol1(a, n)
Invert a vector .
subroutine, public cdiv(a, c, n)
Division of constant c by elements of a .
real(kind=rp), parameter, public neko_m_ln2
subroutine, public chsign(a, n)
Change sign of vector .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public izero(a, n)
Zero an integer vector.
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
subroutine, public copy(a, b, n)
Copy a vector .
subroutine pwmax_scal3(a, b, c, n)
Point-wise maximum of scalar and vector .
subroutine, public add4(a, b, c, d, n)
Vector addition .
pure logical function dabscmp(x, y, tol)
Return double precision absolute comparison .
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine pwmin_vec2(a, b, n)
Point-wise minimum of two vectors .
subroutine pwmin_vec3(a, b, c, n)
Point-wise minimum of two vectors .
subroutine pwmax_vec2(a, b, n)
Point-wise maximum of two vectors .
pure logical function drelcmp(x, y, eps)
Return double precision relative comparison .
real(kind=rp), parameter, public neko_eps
Machine epsilon .
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
pure logical function qrelcmp(x, y, eps)
Return quad precision relative comparison .
subroutine, public rzero(a, n)
Zero a real vector.
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
integer function, public glimax(a, n)
Max of an integer vector of length n.
subroutine sortrp(a, ind, n)
Heap Sort for double precision arrays.
subroutine, public sub2(a, b, n)
Vector substraction .
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
pure logical function srelcmp(x, y, eps)
Return single precision relative comparison .
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
subroutine reordi4(b, ind, n)
reorder single integer array - inverse of swap
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
integer, parameter, public qp
integer, parameter, public i4
integer, parameter, public xp
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.