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)
103 vcross,
vdot2,
vdot3,
vlsc3,
vlsc2,
add2,
add3,
add4,
sub2,
sub3, &
118 real(kind=
sp),
intent(in) :: x
119 real(kind=
sp),
intent(in) :: y
120 real(kind=
sp),
intent(in),
optional :: tol
123 if (
present(tol))
then
133 real(kind=
dp),
intent(in) :: x
134 real(kind=
dp),
intent(in) :: y
135 real(kind=
dp),
intent(in),
optional :: tol
138 if (
present(tol))
then
148 real(kind=
qp),
intent(in) :: x
149 real(kind=
qp),
intent(in) :: y
150 real(kind=
qp),
intent(in),
optional :: tol
153 if (
present(tol))
then
163 real(kind=
sp),
intent(in) :: x
164 real(kind=
sp),
intent(in) :: y
165 real(kind=
sp),
intent(in),
optional :: eps
167 if (
present(eps))
then
168 srelcmp = abs(x - y) .le. eps*abs(y)
177 real(kind=
dp),
intent(in) :: x
178 real(kind=
dp),
intent(in) :: y
179 real(kind=
dp),
intent(in),
optional :: eps
181 if (
present(eps))
then
182 drelcmp = abs(x - y) .le. eps*abs(y)
192 real(kind=
qp),
intent(in) :: x
193 real(kind=
qp),
intent(in) :: y
194 real(kind=
qp),
intent(in),
optional :: eps
196 if (
present(eps))
then
197 qrelcmp = abs(x - y)/abs(y) .lt. eps
206 integer,
intent(in) :: n
207 real(kind=
rp),
dimension(n),
intent(inout) :: a
210 do concurrent(i = 1:n)
217 integer,
intent(in) :: n
218 integer,
dimension(n),
intent(inout) :: a
221 do concurrent(i = 1:n)
228 integer,
intent(in) :: m, n, e
229 real(kind=
rp),
intent(inout) :: a(m,n)
232 do concurrent(j = 1:n)
239 integer,
intent(in) :: n
240 real(kind=
rp),
dimension(n),
intent(inout) :: a
243 do concurrent(i = 1:n)
250 integer,
intent(in) :: n
251 real(kind=
rp),
dimension(n),
intent(in) :: b
252 real(kind=
rp),
dimension(n),
intent(inout) :: a
255 do concurrent(i = 1:n)
269 integer,
intent(in) :: n, n_mask
270 real(kind=
rp),
dimension(n),
intent(in) :: b
271 real(kind=
rp),
dimension(n),
intent(inout) :: a
272 integer,
dimension(0:n_mask) ::
mask
275 do concurrent(i = 1:n_mask)
290 integer,
intent(in) :: n, n_mask
291 real(kind=
rp),
dimension(n),
intent(in) :: b
292 real(kind=
rp),
dimension(n),
intent(inout) :: a
293 integer,
dimension(n_mask) ::
mask
296 do concurrent(i = 1:n_mask)
313 integer,
intent(in) :: n, n_mask
314 real(kind=
rp),
dimension(n),
intent(in) :: b
315 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
316 integer,
dimension(0:n_mask) ::
mask
319 do concurrent(i = 1:n_mask)
336 integer,
intent(in) :: n, n_mask
337 real(kind=
rp),
dimension(n),
intent(in) :: b
338 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
339 integer,
dimension(n_mask) ::
mask
342 do concurrent(i = 1:n_mask)
359 integer,
intent(in) :: n, n_mask
360 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
361 real(kind=
rp),
dimension(n),
intent(inout) :: a
362 integer,
dimension(0:n_mask) ::
mask
365 do concurrent(i = 1:n_mask)
382 integer,
intent(in) :: n, n_mask
383 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
384 real(kind=
rp),
dimension(n),
intent(inout) :: a
385 integer,
dimension(n_mask) ::
mask
388 do concurrent(i = 1:n_mask)
398 integer,
intent(in) :: n, n_mask
399 real(kind=
rp),
dimension(n),
intent(inout) :: a
400 real(kind=
rp),
intent(in) :: c
401 integer,
dimension(n_mask),
intent(in) ::
mask
404 do concurrent(i = 1:n_mask)
412 integer,
intent(in) :: n
413 real(kind=
rp),
dimension(n),
intent(inout) :: a
414 real(kind=
rp),
intent(in) :: c
417 do concurrent(i = 1:n)
424 integer,
intent(in) :: n
425 real(kind=
rp),
dimension(n),
intent(inout) :: a
426 real(kind=
rp),
dimension(n),
intent(in) :: b
427 real(kind=
rp),
intent(in) :: c
430 do concurrent(i = 1:n)
438 integer,
intent(in) :: n
439 real(kind=
rp),
dimension(n),
intent(inout) :: a
440 real(kind=
rp),
intent(in) :: c
443 do concurrent(i = 1:n)
450 integer,
intent(in) :: n
451 real(kind=
rp),
dimension(n),
intent(inout) :: a
452 real(kind=
rp),
dimension(n),
intent(in) :: b
453 real(kind=
rp),
intent(in) :: c
456 do concurrent(i = 1:n)
463 integer,
intent(in) :: n
464 real(kind=
rp),
dimension(n),
intent(inout) :: a
465 real(kind=
rp),
intent(in) :: s
468 do concurrent(i = 1:n)
475 integer,
intent(in) :: n
476 real(kind=
rp),
dimension(n),
intent(inout) :: a
477 real(kind=
rp),
dimension(n),
intent(in) :: b
478 real(kind=
rp),
intent(in) :: s
481 do concurrent(i = 1:n)
488 integer,
intent(in) :: n
489 real(kind=
rp),
dimension(n),
intent(inout) :: a
490 real(kind=
rp),
intent(in) :: c
493 do concurrent(i = 1:n)
500 integer,
intent(in) :: n
501 real(kind=
rp),
dimension(n) :: a
509 call mpi_allreduce(mpi_in_place, tmp, 1, &
517 integer,
intent(in) :: n
518 real(kind=
rp),
dimension(n) :: a
526 call mpi_allreduce(tmp,
glmax, 1, &
532 integer,
intent(in) :: n
533 integer,
dimension(n) :: a
541 call mpi_allreduce(tmp,
glimax, 1, &
547 integer,
intent(in) :: n
548 real(kind=
rp),
dimension(n) :: a
556 call mpi_allreduce(tmp,
glmin, 1, &
562 integer,
intent(in) :: n
563 integer,
dimension(n) :: a
571 call mpi_allreduce(tmp,
glimin, 1, &
580 integer,
intent(in) :: n
581 real(kind=
rp),
dimension(n),
intent(inout) :: a
584 do concurrent(i = 1:n)
593 real(kind=
rp),
intent(in) :: vec(n)
594 real(kind=
rp) :: tmax
597 tmax =
max(tmax, vec(i))
603 integer,
intent(in) :: n
604 real(kind=
rp),
intent(in) :: vec(n)
605 real(kind=
rp) :: tmin
609 tmin = min(tmin, vec(i))
615 integer,
intent(in) :: n
616 real(kind=
rp),
dimension(n),
intent(inout) :: a
619 do concurrent(i = 1:n)
620 a(i) = 1.0_xp /
real(a(i),
xp)
627 integer,
intent(in) :: n
628 real(kind=
rp),
dimension(n),
intent(inout) :: a
629 real(kind=
rp),
dimension(n),
intent(in) :: b, c
632 do concurrent(i = 1:n)
633 a(i) =
real(b(i),
xp) / c(i)
640 integer,
intent(in) :: n
641 real(kind=
rp),
dimension(n),
intent(inout) :: a
642 real(kind=
rp),
dimension(n),
intent(in) :: b
645 do concurrent(i = 1:n)
646 a(i) = 1.0_xp /
real(b(i),
xp)
653 subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
654 integer,
intent(in) :: n
655 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
656 real(kind=
rp),
dimension(n),
intent(in) :: w1, w2, w3
657 real(kind=
rp),
dimension(n),
intent(out) :: u1, u2, u3
660 do concurrent(i = 1:n)
661 u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
662 u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
663 u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
670 subroutine vdot2(dot, u1, u2, v1, v2, n)
671 integer,
intent(in) :: n
672 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2
673 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2
674 real(kind=
rp),
dimension(n),
intent(out) :: dot
676 do concurrent(i = 1:n)
677 dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
684 subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
685 integer,
intent(in) :: n
686 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2, u3
687 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
688 real(kind=
rp),
dimension(n),
intent(out) :: dot
691 do concurrent(i = 1:n)
692 dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
698 function vlsc3(u, v, w, n)
result(s)
699 integer,
intent(in) :: n
700 real(kind=
rp),
dimension(n),
intent(in) :: u, v, w
706 s = s + u(i)*v(i)*w(i)
713 integer,
intent(in) :: n
714 real(kind=
rp),
dimension(n),
intent(in) :: u, v
727 integer,
intent(in) :: n
728 real(kind=
rp),
dimension(n),
intent(inout) :: a
729 real(kind=
rp),
dimension(n),
intent(in) :: b
732 do concurrent(i = 1:n)
740 integer,
intent(in) :: n
741 real(kind=
rp),
dimension(n),
intent(inout) :: a
742 real(kind=
rp),
dimension(n),
intent(in) :: b
743 real(kind=
rp),
dimension(n),
intent(in) :: c
746 do concurrent(i = 1:n)
754 integer,
intent(in) :: n
755 real(kind=
rp),
dimension(n),
intent(out) :: a
756 real(kind=
rp),
dimension(n),
intent(in) :: d
757 real(kind=
rp),
dimension(n),
intent(in) :: c
758 real(kind=
rp),
dimension(n),
intent(in) :: b
761 do concurrent(i = 1:n)
762 a(i) = b(i) + c(i) + d(i)
769 integer,
intent(in) :: n
770 real(kind=
rp),
dimension(n),
intent(inout) :: a
771 real(kind=
rp),
dimension(n),
intent(in) :: b
774 do concurrent(i = 1:n)
782 integer,
intent(in) :: n
783 real(kind=
rp),
dimension(n),
intent(inout) :: a
784 real(kind=
rp),
dimension(n),
intent(in) :: b
785 real(kind=
rp),
dimension(n),
intent(in) :: c
788 do concurrent(i = 1:n)
798 integer,
intent(in) :: n
799 real(kind=
rp),
dimension(n),
intent(inout) :: a
800 real(kind=
rp),
dimension(n),
intent(in) :: b
801 real(kind=
rp),
intent(in) :: c1
804 do concurrent(i = 1:n)
805 a(i) = c1 * a(i) + b(i)
813 integer,
intent(in) :: n
814 real(kind=
rp),
dimension(n),
intent(inout) :: a
815 real(kind=
rp),
dimension(n),
intent(in) :: b
816 real(kind=
rp),
intent(in) :: c1
819 do concurrent(i = 1:n)
820 a(i) = a(i) + c1 * b(i)
827 integer,
intent(in) :: n
828 real(kind=
rp),
dimension(n),
intent(inout) :: a
829 real(kind=
rp),
dimension(n),
intent(in) :: b
830 real(kind=
rp),
intent(in) :: c1
833 do concurrent(i = 1:n)
834 a(i) = a(i) + c1 * ( b(i) * b(i) )
841 integer,
intent(in) :: n
842 real(kind=
rp),
dimension(n),
intent(inout) :: a
843 real(kind=
rp),
dimension(n),
intent(in) :: b
846 do concurrent(i = 1:n)
847 a(i) =
real(a(i),
xp) / b(i)
855 integer,
intent(in) :: n
856 real(kind=
rp),
dimension(n),
intent(inout) :: a
857 real(kind=
rp),
dimension(n),
intent(in) :: b
860 do concurrent(i = 1:n)
868 integer,
intent(in) :: n
869 real(kind=
rp),
dimension(n),
intent(inout) :: a
870 real(kind=
rp),
dimension(n),
intent(in) :: b
871 real(kind=
rp),
dimension(n),
intent(in) :: c
874 do concurrent(i = 1:n)
882 integer,
intent(in) :: n
883 real(kind=
rp),
dimension(n),
intent(inout) :: a
884 real(kind=
rp),
dimension(n),
intent(in) :: b
885 real(kind=
rp),
dimension(n),
intent(in) :: c
888 do concurrent(i = 1:n)
889 a(i) = a(i) - b(i) * c(i)
896 integer,
intent(in) :: n
897 real(kind=
rp),
dimension(n),
intent(inout) :: a
898 real(kind=
rp),
dimension(n),
intent(in) :: b
899 real(kind=
rp),
dimension(n),
intent(in) :: c
900 real(kind=
rp),
intent(in) :: c1, c2
903 do concurrent(i = 1:n)
904 a(i) = c1 * b(i) + c2 * c(i)
910 subroutine add4s3(a, b, c, d, c1, c2, c3, n)
911 integer,
intent(in) :: n
912 real(kind=
rp),
dimension(n),
intent(inout) :: a
913 real(kind=
rp),
dimension(n),
intent(in) :: b
914 real(kind=
rp),
dimension(n),
intent(in) :: c
915 real(kind=
rp),
dimension(n),
intent(in) :: d
916 real(kind=
rp),
intent(in) :: c1, c2, c3
919 do concurrent(i = 1:n)
920 a(i) = c1 * b(i) + c2 * c(i) + c3 * d(i)
926 subroutine add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
927 integer,
intent(in) :: n
928 real(kind=
rp),
dimension(n),
intent(inout) :: a
929 real(kind=
rp),
dimension(n),
intent(in) :: b
930 real(kind=
rp),
dimension(n),
intent(in) :: c
931 real(kind=
rp),
dimension(n),
intent(in) :: d
932 real(kind=
rp),
dimension(n),
intent(in) :: e
933 real(kind=
rp),
intent(in) :: c1, c2, c3, c4
936 do concurrent(i = 1:n)
937 a(i) = a(i) + c1 * b(i) + c2 * c(i) + c3 * d(i) + c4 * e(i)
944 integer,
intent(in) :: n
945 real(kind=
rp),
dimension(n),
intent(inout) :: a
946 real(kind=
rp),
dimension(n),
intent(in) :: b
947 real(kind=
rp),
dimension(n),
intent(in) :: c
948 real(kind=
rp),
dimension(n),
intent(in) :: d
951 do concurrent(i = 1:n)
952 a(i) = a(i) - b(i) * c(i) * d(i)
959 integer,
intent(in) :: n
960 real(kind=
rp),
dimension(n),
intent(inout) :: a
961 real(kind=
rp),
dimension(n),
intent(in) :: b
962 real(kind=
rp),
dimension(n),
intent(in) :: c
965 do concurrent(i = 1:n)
966 a(i) = a(i) + b(i) * c(i)
973 integer,
intent(in) :: n
974 real(kind=
rp),
dimension(n),
intent(inout) :: a
975 real(kind=
rp),
dimension(n),
intent(in) :: b
976 real(kind=
rp),
dimension(n),
intent(in) :: c
977 real(kind=
rp),
dimension(n),
intent(in) :: d
980 do concurrent(i = 1:n)
981 a(i) = a(i) + b(i) * c(i) * d(i)
988 integer,
intent(in) :: n
989 real(kind=
rp),
dimension(n),
intent(inout) :: a
990 real(kind=
rp),
dimension(n),
intent(in) :: b
991 real(kind=
rp),
dimension(n),
intent(in) :: c
992 real(kind=
rp),
intent(in) :: s
995 do concurrent(i = 1:n)
996 a(i) = a(i) + s * b(i) * c(i)
1003 integer,
intent(in) :: n
1004 real(kind=
rp),
dimension(n),
intent(inout) :: a
1005 real(kind=
rp),
dimension(n),
intent(in) :: b
1006 real(kind=
rp),
dimension(n),
intent(in) :: c
1007 real(kind=
rp),
dimension(n),
intent(in) :: d
1008 real(kind=
rp),
dimension(n),
intent(in) :: e
1011 do concurrent(i = 1:n)
1012 a(i) = b(i)*c(i) - d(i)*e(i)
1019 integer,
intent(in) :: n
1020 real(kind=
rp),
dimension(n),
intent(inout) :: a
1021 real(kind=
rp),
dimension(n),
intent(in) :: b
1022 real(kind=
rp),
dimension(n),
intent(in) :: c
1023 real(kind=
rp),
intent(in) :: c1, c2
1026 do concurrent(i = 1:n)
1027 a(i) = b(i) + c1*(a(i)-c2*c(i))
1034 integer,
intent(in) :: n
1035 real(kind=
rp),
dimension(n),
intent(inout) :: a
1036 real(kind=
rp),
dimension(n),
intent(in) :: b
1037 real(kind=
rp),
dimension(n),
intent(in) :: c
1038 real(kind=
rp),
intent(in) :: c1, c2
1041 do concurrent(i = 1:n)
1042 a(i) = a(i) + c1*b(i)+c2*c(i)
1049 integer,
intent(in) :: n
1050 real(kind=
rp),
dimension(n),
intent(in) :: a
1051 real(kind=
rp),
dimension(n),
intent(in) :: b
1053 real(kind=
xp) :: tmp
1058 tmp = tmp + a(i) * b(i)
1061 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
1071 real(kind=
rp),
dimension(n),
intent(in) :: c
1073 real(kind=
xp) :: tmp
1078 tmp = tmp + a(i) * b(i) * c(i)
1081 call mpi_allreduce(mpi_in_place, tmp, 1, &
1087 integer,
intent(in) :: n
1088 real(kind=
rp),
dimension(n),
intent(in) :: a
1089 real(kind=
rp),
dimension(n),
intent(in) :: b
1090 real(kind=
rp),
dimension(n),
intent(in) :: c
1091 real(kind=
rp),
dimension(n),
intent(in) :: d
1093 real(kind=
xp) :: tmp
1098 tmp = tmp + a(i) * b(i) * c(i) * d(i)
1101 call mpi_allreduce(mpi_in_place, tmp, 1, &
1110 integer,
intent(in) :: n
1111 real(kind=
rp),
dimension(n),
intent(in) :: a
1112 real(kind=
rp),
dimension(n),
intent(in) :: b
1114 real(kind=
xp) :: tmp
1119 tmp = tmp + (a(i) - b(i))**2
1122 call mpi_allreduce(mpi_in_place, tmp, 1, &
1134 integer,
intent(in) :: n
1135 real(kind=
rp),
intent(inout) :: a(n)
1136 integer,
intent(out) :: ind(n)
1138 integer :: j, ir, i, ii, l
1144 if (n .le. 1)
return
1168 do while (j .le. ir)
1170 if ( a(j) .lt. a(j+1) ) j = j + 1
1172 if (aa .lt. a(j))
then
1192 integer,
intent(in) :: n
1193 integer(i4),
intent(inout) :: a(n)
1194 integer,
intent(out) :: ind(n)
1196 integer :: j, ir, i, ii, l
1202 if (n .le. 1)
return
1225 do while (j .le. ir)
1227 if ( a(j) .lt. a(j + 1) ) j = j + 1
1229 if (aa .lt. a(j))
then
1248 integer,
intent(in) :: n
1249 real(kind=
rp),
intent(inout) :: b(n)
1250 integer,
intent(in) :: ind(n)
1251 real(kind=
rp) :: temp(n)
1268 integer,
intent(in) :: n
1269 integer(i4),
intent(inout) :: b(n)
1270 integer,
intent(in) :: ind(n)
1271 integer(i4) :: temp(n)
1288 integer,
intent(in) :: n
1289 real(kind=
rp),
intent(inout) :: b(n)
1290 integer,
intent(in) :: ind(n)
1291 real(kind=
rp) :: temp(n)
1308 integer,
intent(in) :: n
1309 integer(i4),
intent(inout) :: b(n)
1310 integer,
intent(in) :: ind(n)
1311 integer(i4) :: temp(n)
1328 integer,
intent(in) :: n
1329 real(kind=
rp),
intent(inout) :: b(n)
1330 integer,
intent(inout) :: ind(n)
1331 real(kind=
rp) :: temp(n)
1332 integer :: tempind(n)
1338 tempind(jj) = ind(i)
1351 integer,
intent(in) :: n
1352 integer(i4),
intent(inout) :: b(n)
1353 integer,
intent(inout) :: ind(n)
1354 integer(i4) :: temp(n)
1355 integer :: tempind(n)
1361 tempind(jj) = ind(i)
1373 integer,
intent(in) :: n
1374 real(kind=
rp),
dimension(n),
intent(inout) :: a
1386 integer,
intent(in) :: n
1387 real(kind=
rp),
dimension(n),
intent(inout) :: a
1388 real(kind=
rp),
dimension(n),
intent(in) :: b
1392 a(i) =
max(a(i), b(i))
1398 integer,
intent(in) :: n
1399 real(kind=
rp),
dimension(n),
intent(inout) :: a
1400 real(kind=
rp),
dimension(n),
intent(in) :: b, c
1404 a(i) =
max(b(i), c(i))
1410 integer,
intent(in) :: n
1411 real(kind=
rp),
dimension(n),
intent(inout) :: a
1412 real(kind=
rp),
intent(in) :: b
1422 integer,
intent(in) :: n
1423 real(kind=
rp),
dimension(n),
intent(inout) :: a
1424 real(kind=
rp),
dimension(n),
intent(in) :: b
1425 real(kind=
rp),
intent(in) :: c
1435 integer,
intent(in) :: n
1436 real(kind=
rp),
dimension(n),
intent(inout) :: a
1437 real(kind=
rp),
dimension(n),
intent(in) :: b
1441 a(i) = min(a(i), b(i))
1447 integer,
intent(in) :: n
1448 real(kind=
rp),
dimension(n),
intent(inout) :: a
1449 real(kind=
rp),
dimension(n),
intent(in) :: b, c
1453 a(i) = min(b(i), c(i))
1459 integer,
intent(in) :: n
1460 real(kind=
rp),
dimension(n),
intent(inout) :: a
1461 real(kind=
rp),
intent(in) :: b
1471 integer,
intent(in) :: n
1472 real(kind=
rp),
dimension(n),
intent(inout) :: a
1473 real(kind=
rp),
dimension(n),
intent(in) :: b
1474 real(kind=
rp),
intent(in) :: c
1485 function matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33) &
1487 real(kind=
rp),
intent(in) :: a11, a12, a13, a21, a22, a23, a31, a32, a33
1508 real(kind=
xp),
intent(in) :: a(3,3)
1509 real(kind=
xp) :: b(3,3)
1510 real(kind=
xp) :: detinv
1514 detinv = 1.0_xp /
real(a(1,1)*a(2,2)*a(3,3) - a(1,1)*a(2,3)*a(3,2) &
1515 - a(1,2)*a(2,1)*a(3,3) + a(1,2)*a(2,3)*a(3,1)&
1516 + a(1,3)*a(2,1)*a(3,2) - a(1,3)*a(2,2)*a(3,1),
xp)
1519 b(1,1) = +detinv * (a(2,2)*a(3,3) - a(2,3)*a(3,2))
1520 b(2,1) = -detinv * (a(2,1)*a(3,3) - a(2,3)*a(3,1))
1521 b(3,1) = +detinv * (a(2,1)*a(3,2) - a(2,2)*a(3,1))
1522 b(1,2) = -detinv * (a(1,2)*a(3,3) - a(1,3)*a(3,2))
1523 b(2,2) = +detinv * (a(1,1)*a(3,3) - a(1,3)*a(3,1))
1524 b(3,2) = -detinv * (a(1,1)*a(3,2) - a(1,2)*a(3,1))
1525 b(1,3) = +detinv * (a(1,2)*a(2,3) - a(1,3)*a(2,2))
1526 b(2,3) = -detinv * (a(1,1)*a(2,3) - a(1,3)*a(2,1))
1527 b(3,3) = +detinv * (a(1,1)*a(2,2) - a(1,2)*a(2,1))
1533 real(kind=
rp),
intent(in) :: x
1534 real(kind=
rp) :: val
1535 real(kind=
rp),
parameter :: xdmin = 0.0001_rp
1536 real(kind=
rp),
parameter :: xdmax = 0.9999_rp
1539 if (x <= xdmin)
then
1542 else if (x >= xdmax)
then
1547 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1550 val = 1.0_rp / (1.0_rp + exp(g))
1556 real(kind=
rp),
intent(in) :: x
1557 real(kind=
rp) :: val
1558 real(kind=
rp),
parameter :: xdmin = 0.0001_rp
1559 real(kind=
rp),
parameter :: xdmax = 0.9999_rp
1560 real(kind=
rp) :: arg, g, dg, s_val
1562 if (x <= xdmin .or. x >= xdmax)
then
1569 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1572 dg = -(1.0_rp / ((x - 1.0_rp)**2)) - (1.0_rp / (x**2))
1575 s_val = 1.0_rp / (1.0_rp + exp(g))
1577 val = -s_val * (1.0_rp - s_val) * dg
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, public dabscmp(x, y, tol)
Return double precision absolute comparison .
real(kind=rp), parameter, public pi
pure logical function qabscmp(x, y, tol)
Return double precision absolute comparison .
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 addcol3s2(a, b, c, s, 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 .
real(rp) function, dimension(3, 3), public matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33)
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 .
real(kind=rp) function, public math_stepf(x)
Smooth step function S(x) Returns 0 for x <= 0, 1 for x >= 1, and smooth transition in between.
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)
subroutine, public cpwmin2(a, b, n)
Point-wise minimum of scalar and vector .
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 cpwmin3(a, b, c, n)
Point-wise minimum of scalar and vector .
subroutine, public pwmax3(a, b, c, n)
Point-wise maximum of two vectors .
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 .
real(kind=xp) function, dimension(3, 3), public matinv3(a)
Performs a direct calculation of the inverse of a 3×3 matrix. M33INV and M44INV by David G....
subroutine, public pwmax2(a, b, n)
Point-wise maximum of two vectors .
subroutine, public pwmin2(a, b, n)
Point-wise minimum of two vectors .
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 cpwmax3(a, b, c, n)
Point-wise maximum of scalar and 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, public add4s3(a, b, c, d, c1, c2, c3, n)
Returns .
subroutine, public add4(a, b, c, d, n)
Vector addition .
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
Returns .
real(kind=rp) function, public math_dstepf(x)
Derivative of math_stepf with respect to x: d(stepf)/dx.
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, public sabscmp(x, y, tol)
Return single precision absolute comparison .
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.
subroutine, public cpwmax2(a, b, n)
Point-wise maximum of scalar and vector .
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 pwmin3(a, b, c, n)
Point-wise minimum of two vectors .
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.