63 use mpi_f08,
only : mpi_min, mpi_max, mpi_sum, mpi_in_place, mpi_integer, &
70 real(kind=
rp),
public,
parameter ::
neko_eps = epsilon(1.0_rp)
76 real(kind=
rp),
public,
parameter ::
pi = 4._rp*atan(1._rp)
104 vcross,
vdot2,
vdot3,
vlsc3,
vlsc2,
add2,
add3,
add4,
sub2,
sub3, &
120 real(kind=
sp),
intent(in) :: x
121 real(kind=
sp),
intent(in) :: y
122 real(kind=
sp),
intent(in),
optional :: tol
125 if (
present(tol))
then
135 real(kind=
dp),
intent(in) :: x
136 real(kind=
dp),
intent(in) :: y
137 real(kind=
dp),
intent(in),
optional :: tol
140 if (
present(tol))
then
150 real(kind=
qp),
intent(in) :: x
151 real(kind=
qp),
intent(in) :: y
152 real(kind=
qp),
intent(in),
optional :: tol
155 if (
present(tol))
then
165 real(kind=
sp),
intent(in) :: x
166 real(kind=
sp),
intent(in) :: y
167 real(kind=
sp),
intent(in),
optional :: eps
169 if (
present(eps))
then
170 srelcmp = abs(x - y) .le. eps*abs(y)
179 real(kind=
dp),
intent(in) :: x
180 real(kind=
dp),
intent(in) :: y
181 real(kind=
dp),
intent(in),
optional :: eps
183 if (
present(eps))
then
184 drelcmp = abs(x - y) .le. eps*abs(y)
194 real(kind=
qp),
intent(in) :: x
195 real(kind=
qp),
intent(in) :: y
196 real(kind=
qp),
intent(in),
optional :: eps
198 if (
present(eps))
then
199 qrelcmp = abs(x - y)/abs(y) .lt. eps
208 integer,
intent(in) :: n
209 real(kind=
rp),
dimension(n),
intent(inout) :: a
212 do concurrent(i = 1:n)
219 integer,
intent(in) :: n
220 integer,
dimension(n),
intent(inout) :: a
223 do concurrent(i = 1:n)
230 integer,
intent(in) :: m, n, e
231 real(kind=
rp),
intent(inout) :: a(m,n)
234 do concurrent(j = 1:n)
241 integer,
intent(in) :: n
242 real(kind=
rp),
dimension(n),
intent(inout) :: a
245 do concurrent(i = 1:n)
252 integer,
intent(in) :: n
253 real(kind=
rp),
dimension(n),
intent(in) :: b
254 real(kind=
rp),
dimension(n),
intent(inout) :: a
257 do concurrent(i = 1:n)
271 integer,
intent(in) :: n, n_mask
272 real(kind=
rp),
dimension(n),
intent(in) :: b
273 real(kind=
rp),
dimension(n),
intent(inout) :: a
274 integer,
dimension(0:n_mask) ::
mask
277 do concurrent(i = 1:n_mask)
292 integer,
intent(in) :: n, n_mask
293 real(kind=
rp),
dimension(n),
intent(in) :: b
294 real(kind=
rp),
dimension(n),
intent(inout) :: a
295 integer,
dimension(n_mask) ::
mask
298 do concurrent(i = 1:n_mask)
315 integer,
intent(in) :: n, n_mask
316 real(kind=
rp),
dimension(n),
intent(in) :: b
317 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
318 integer,
dimension(0:n_mask) ::
mask
321 do concurrent(i = 1:n_mask)
338 integer,
intent(in) :: lx, ly, lz, n_mask
339 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
340 real(kind=
rp),
dimension(:, :, :, :),
intent(in) :: b
341 integer,
dimension(0:n_mask),
intent(in) ::
mask
342 integer,
dimension(0:n_mask),
intent(in) :: facet
349 select case (facet(l))
351 a(l) = b(idx(2), idx(3), facet(l), idx(4))
353 a(l) = b(idx(1), idx(3), facet(l), idx(4))
355 a(l) = b(idx(1), idx(2), facet(l), idx(4))
371 integer,
intent(in) :: n, n_mask
372 real(kind=
rp),
dimension(n),
intent(in) :: b
373 real(kind=
rp),
dimension(n_mask),
intent(inout) :: a
374 integer,
dimension(n_mask) ::
mask
377 do concurrent(i = 1:n_mask)
394 integer,
intent(in) :: n, n_mask
395 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
396 real(kind=
rp),
dimension(n),
intent(inout) :: a
397 integer,
dimension(0:n_mask) ::
mask
400 do concurrent(i = 1:n_mask)
417 integer,
intent(in) :: n, n_mask
418 real(kind=
rp),
dimension(n_mask),
intent(in) :: b
419 real(kind=
rp),
dimension(n),
intent(inout) :: a
420 integer,
dimension(n_mask) ::
mask
423 do concurrent(i = 1:n_mask)
433 integer,
intent(in) :: n, n_mask
434 real(kind=
rp),
dimension(n),
intent(inout) :: a
435 real(kind=
rp),
intent(in) :: c
436 integer,
dimension(n_mask),
intent(in) ::
mask
439 do concurrent(i = 1:n_mask)
447 integer,
intent(in) :: n
448 real(kind=
rp),
dimension(n),
intent(inout) :: a
449 real(kind=
rp),
intent(in) :: c
452 do concurrent(i = 1:n)
459 integer,
intent(in) :: n
460 real(kind=
rp),
dimension(n),
intent(inout) :: a
461 real(kind=
rp),
dimension(n),
intent(in) :: b
462 real(kind=
rp),
intent(in) :: c
465 do concurrent(i = 1:n)
473 integer,
intent(in) :: n
474 real(kind=
rp),
dimension(n),
intent(inout) :: a
475 real(kind=
rp),
intent(in) :: c
478 do concurrent(i = 1:n)
485 integer,
intent(in) :: n
486 real(kind=
rp),
dimension(n),
intent(inout) :: a
487 real(kind=
rp),
dimension(n),
intent(in) :: b
488 real(kind=
rp),
intent(in) :: c
491 do concurrent(i = 1:n)
498 integer,
intent(in) :: n
499 real(kind=
rp),
dimension(n),
intent(inout) :: a
500 real(kind=
rp),
intent(in) :: s
503 do concurrent(i = 1:n)
510 integer,
intent(in) :: n
511 real(kind=
rp),
dimension(n),
intent(inout) :: a
512 real(kind=
rp),
dimension(n),
intent(in) :: b
513 real(kind=
rp),
intent(in) :: s
516 do concurrent(i = 1:n)
523 integer,
intent(in) :: n
524 real(kind=
rp),
dimension(n),
intent(inout) :: a
525 real(kind=
rp),
intent(in) :: c
528 do concurrent(i = 1:n)
534 subroutine cwrap(a, min_val, max_val, n)
535 integer,
intent(in) :: n
536 real(kind=
rp),
dimension(n),
intent(inout) :: a
537 real(kind=
rp),
intent(in) :: min_val, max_val
540 if (n .lt. 1 .or. max_val .le. min_val)
return
542 do concurrent(i = 1:n)
543 a(i) = modulo(a(i) - min_val, max_val - min_val) + min_val
549 integer,
intent(in) :: n
550 real(kind=
rp),
dimension(n) :: a
558 call mpi_allreduce(mpi_in_place, tmp, 1, &
566 integer,
intent(in) :: n
567 real(kind=
rp),
dimension(n) :: a
575 call mpi_allreduce(tmp,
glmax, 1, &
581 integer,
intent(in) :: n
582 integer,
dimension(n) :: a
590 call mpi_allreduce(tmp,
glimax, 1, &
596 integer,
intent(in) :: n
597 real(kind=
rp),
dimension(n) :: a
605 call mpi_allreduce(tmp,
glmin, 1, &
611 integer,
intent(in) :: n
612 integer,
dimension(n) :: a
620 call mpi_allreduce(tmp,
glimin, 1, &
629 integer,
intent(in) :: n
630 real(kind=
rp),
dimension(n),
intent(inout) :: a
633 do concurrent(i = 1:n)
642 real(kind=
rp),
intent(in) :: vec(n)
643 real(kind=
rp) :: tmax
646 tmax =
max(tmax, vec(i))
652 integer,
intent(in) :: n
653 real(kind=
rp),
intent(in) :: vec(n)
654 real(kind=
rp) :: tmin
658 tmin = min(tmin, vec(i))
664 integer,
intent(in) :: n
665 real(kind=
rp),
dimension(n),
intent(inout) :: a
668 do concurrent(i = 1:n)
669 a(i) = 1.0_xp /
real(a(i),
xp)
676 integer,
intent(in) :: n
677 real(kind=
rp),
dimension(n),
intent(inout) :: a
678 real(kind=
rp),
dimension(n),
intent(in) :: b, c
681 do concurrent(i = 1:n)
682 a(i) =
real(b(i),
xp) / c(i)
689 integer,
intent(in) :: n
690 real(kind=
rp),
dimension(n),
intent(inout) :: a
691 real(kind=
rp),
dimension(n),
intent(in) :: b
694 do concurrent(i = 1:n)
695 a(i) = 1.0_xp /
real(b(i),
xp)
702 subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
703 integer,
intent(in) :: n
704 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
705 real(kind=
rp),
dimension(n),
intent(in) :: w1, w2, w3
706 real(kind=
rp),
dimension(n),
intent(out) :: u1, u2, u3
709 do concurrent(i = 1:n)
710 u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
711 u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
712 u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
719 subroutine vdot2(dot, u1, u2, v1, v2, n)
720 integer,
intent(in) :: n
721 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2
722 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2
723 real(kind=
rp),
dimension(n),
intent(out) :: dot
725 do concurrent(i = 1:n)
726 dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
733 subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
734 integer,
intent(in) :: n
735 real(kind=
rp),
dimension(n),
intent(in) :: u1, u2, u3
736 real(kind=
rp),
dimension(n),
intent(in) :: v1, v2, v3
737 real(kind=
rp),
dimension(n),
intent(out) :: dot
740 do concurrent(i = 1:n)
741 dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
747 function vlsc3(u, v, w, n)
result(s)
748 integer,
intent(in) :: n
749 real(kind=
rp),
dimension(n),
intent(in) :: u, v, w
755 s = s + u(i)*v(i)*w(i)
762 integer,
intent(in) :: n
763 real(kind=
rp),
dimension(n),
intent(in) :: u, v
776 integer,
intent(in) :: n
777 real(kind=
rp),
dimension(n),
intent(inout) :: a
778 real(kind=
rp),
dimension(n),
intent(in) :: b
781 do concurrent(i = 1:n)
789 integer,
intent(in) :: n
790 real(kind=
rp),
dimension(n),
intent(inout) :: a
791 real(kind=
rp),
dimension(n),
intent(in) :: b
792 real(kind=
rp),
dimension(n),
intent(in) :: c
795 do concurrent(i = 1:n)
803 integer,
intent(in) :: n
804 real(kind=
rp),
dimension(n),
intent(out) :: a
805 real(kind=
rp),
dimension(n),
intent(in) :: d
806 real(kind=
rp),
dimension(n),
intent(in) :: c
807 real(kind=
rp),
dimension(n),
intent(in) :: b
810 do concurrent(i = 1:n)
811 a(i) = b(i) + c(i) + d(i)
818 integer,
intent(in) :: n
819 real(kind=
rp),
dimension(n),
intent(inout) :: a
820 real(kind=
rp),
dimension(n),
intent(in) :: b
823 do concurrent(i = 1:n)
831 integer,
intent(in) :: n
832 real(kind=
rp),
dimension(n),
intent(inout) :: a
833 real(kind=
rp),
dimension(n),
intent(in) :: b
834 real(kind=
rp),
dimension(n),
intent(in) :: c
837 do concurrent(i = 1:n)
847 integer,
intent(in) :: n
848 real(kind=
rp),
dimension(n),
intent(inout) :: a
849 real(kind=
rp),
dimension(n),
intent(in) :: b
850 real(kind=
rp),
intent(in) :: c1
853 do concurrent(i = 1:n)
854 a(i) = c1 * a(i) + b(i)
862 integer,
intent(in) :: n
863 real(kind=
rp),
dimension(n),
intent(inout) :: a
864 real(kind=
rp),
dimension(n),
intent(in) :: b
865 real(kind=
rp),
intent(in) :: c1
868 do concurrent(i = 1:n)
869 a(i) = a(i) + c1 * b(i)
876 integer,
intent(in) :: n
877 real(kind=
rp),
dimension(n),
intent(inout) :: a
878 real(kind=
rp),
dimension(n),
intent(in) :: b
879 real(kind=
rp),
intent(in) :: c1
882 do concurrent(i = 1:n)
883 a(i) = a(i) + c1 * ( b(i) * b(i) )
890 integer,
intent(in) :: n
891 real(kind=
rp),
dimension(n),
intent(inout) :: a
892 real(kind=
rp),
dimension(n),
intent(in) :: b
895 do concurrent(i = 1:n)
896 a(i) =
real(a(i),
xp) / b(i)
904 integer,
intent(in) :: n
905 real(kind=
rp),
dimension(n),
intent(inout) :: a
906 real(kind=
rp),
dimension(n),
intent(in) :: b
909 do concurrent(i = 1:n)
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
923 do concurrent(i = 1:n)
931 integer,
intent(in) :: n
932 real(kind=
rp),
dimension(n),
intent(inout) :: a
933 real(kind=
rp),
dimension(n),
intent(in) :: b
934 real(kind=
rp),
dimension(n),
intent(in) :: c
937 do concurrent(i = 1:n)
938 a(i) = a(i) - b(i) * c(i)
945 integer,
intent(in) :: n
946 real(kind=
rp),
dimension(n),
intent(inout) :: a
947 real(kind=
rp),
dimension(n),
intent(in) :: b
948 real(kind=
rp),
dimension(n),
intent(in) :: c
949 real(kind=
rp),
intent(in) :: c1, c2
952 do concurrent(i = 1:n)
953 a(i) = c1 * b(i) + c2 * c(i)
959 subroutine add4s3(a, b, c, d, c1, c2, c3, n)
960 integer,
intent(in) :: n
961 real(kind=
rp),
dimension(n),
intent(inout) :: a
962 real(kind=
rp),
dimension(n),
intent(in) :: b
963 real(kind=
rp),
dimension(n),
intent(in) :: c
964 real(kind=
rp),
dimension(n),
intent(in) :: d
965 real(kind=
rp),
intent(in) :: c1, c2, c3
968 do concurrent(i = 1:n)
969 a(i) = c1 * b(i) + c2 * c(i) + c3 * d(i)
975 subroutine add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
976 integer,
intent(in) :: n
977 real(kind=
rp),
dimension(n),
intent(inout) :: a
978 real(kind=
rp),
dimension(n),
intent(in) :: b
979 real(kind=
rp),
dimension(n),
intent(in) :: c
980 real(kind=
rp),
dimension(n),
intent(in) :: d
981 real(kind=
rp),
dimension(n),
intent(in) :: e
982 real(kind=
rp),
intent(in) :: c1, c2, c3, c4
985 do concurrent(i = 1:n)
986 a(i) = a(i) + c1 * b(i) + c2 * c(i) + c3 * d(i) + c4 * e(i)
993 integer,
intent(in) :: n
994 real(kind=
rp),
dimension(n),
intent(inout) :: a
995 real(kind=
rp),
dimension(n),
intent(in) :: b
996 real(kind=
rp),
dimension(n),
intent(in) :: c
997 real(kind=
rp),
dimension(n),
intent(in) :: d
1000 do concurrent(i = 1:n)
1001 a(i) = a(i) - b(i) * c(i) * d(i)
1008 integer,
intent(in) :: n
1009 real(kind=
rp),
dimension(n),
intent(inout) :: a
1010 real(kind=
rp),
dimension(n),
intent(in) :: b
1011 real(kind=
rp),
dimension(n),
intent(in) :: c
1014 do concurrent(i = 1:n)
1015 a(i) = a(i) + b(i) * c(i)
1022 integer,
intent(in) :: n
1023 real(kind=
rp),
dimension(n),
intent(inout) :: a
1024 real(kind=
rp),
dimension(n),
intent(in) :: b
1025 real(kind=
rp),
dimension(n),
intent(in) :: c
1026 real(kind=
rp),
dimension(n),
intent(in) :: d
1029 do concurrent(i = 1:n)
1030 a(i) = a(i) + b(i) * c(i) * d(i)
1037 integer,
intent(in) :: n
1038 real(kind=
rp),
dimension(n),
intent(inout) :: a
1039 real(kind=
rp),
dimension(n),
intent(in) :: b
1040 real(kind=
rp),
dimension(n),
intent(in) :: c
1041 real(kind=
rp),
intent(in) :: s
1044 do concurrent(i = 1:n)
1045 a(i) = a(i) + s * b(i) * c(i)
1052 integer,
intent(in) :: n
1053 real(kind=
rp),
dimension(n),
intent(inout) :: a
1054 real(kind=
rp),
dimension(n),
intent(in) :: b
1055 real(kind=
rp),
dimension(n),
intent(in) :: c
1056 real(kind=
rp),
dimension(n),
intent(in) :: d
1057 real(kind=
rp),
dimension(n),
intent(in) :: e
1060 do concurrent(i = 1:n)
1061 a(i) = b(i)*c(i) - d(i)*e(i)
1068 integer,
intent(in) :: n
1069 real(kind=
rp),
dimension(n),
intent(inout) :: a
1070 real(kind=
rp),
dimension(n),
intent(in) :: b
1071 real(kind=
rp),
dimension(n),
intent(in) :: c
1072 real(kind=
rp),
intent(in) :: c1, c2
1075 do concurrent(i = 1:n)
1076 a(i) = b(i) + c1*(a(i)-c2*c(i))
1083 integer,
intent(in) :: n
1084 real(kind=
rp),
dimension(n),
intent(inout) :: a
1085 real(kind=
rp),
dimension(n),
intent(in) :: b
1086 real(kind=
rp),
dimension(n),
intent(in) :: c
1087 real(kind=
rp),
intent(in) :: c1, c2
1090 do concurrent(i = 1:n)
1091 a(i) = a(i) + c1*b(i)+c2*c(i)
1098 integer,
intent(in) :: n
1099 real(kind=
rp),
dimension(n),
intent(in) :: a
1100 real(kind=
rp),
dimension(n),
intent(in) :: b
1102 real(kind=
xp) :: tmp
1107 tmp = tmp + a(i) * b(i)
1110 call mpi_allreduce(mpi_in_place, tmp, 1, &
1117 integer,
intent(in) :: n
1118 real(kind=
rp),
dimension(n),
intent(in) :: a
1119 real(kind=
rp),
dimension(n),
intent(in) :: b
1120 real(kind=
rp),
dimension(n),
intent(in) :: c
1122 real(kind=
xp) :: tmp
1127 tmp = tmp + a(i) * b(i) * c(i)
1130 call mpi_allreduce(mpi_in_place, tmp, 1, &
1136 integer,
intent(in) :: n
1137 real(kind=
rp),
dimension(n),
intent(in) :: a
1138 real(kind=
rp),
dimension(n),
intent(in) :: b
1139 real(kind=
rp),
dimension(n),
intent(in) :: c
1140 real(kind=
rp),
dimension(n),
intent(in) :: d
1142 real(kind=
xp) :: tmp
1147 tmp = tmp + a(i) * b(i) * c(i) * d(i)
1150 call mpi_allreduce(mpi_in_place, tmp, 1, &
1159 integer,
intent(in) :: n
1160 real(kind=
rp),
dimension(n),
intent(in) :: a
1161 real(kind=
rp),
dimension(n),
intent(in) :: b
1163 real(kind=
xp) :: tmp
1168 tmp = tmp + (a(i) - b(i))**2
1171 call mpi_allreduce(mpi_in_place, tmp, 1, &
1183 integer,
intent(in) :: n
1184 real(kind=
rp),
intent(inout) :: a(n)
1185 integer,
intent(out) :: ind(n)
1187 integer :: j, ir, i, ii, l
1193 if (n .le. 1)
return
1217 do while (j .le. ir)
1219 if ( a(j) .lt. a(j+1) ) j = j + 1
1221 if (aa .lt. a(j))
then
1241 integer,
intent(in) :: n
1242 integer(i4),
intent(inout) :: a(n)
1243 integer,
intent(out) :: ind(n)
1245 integer :: j, ir, i, ii, l
1251 if (n .le. 1)
return
1274 do while (j .le. ir)
1276 if ( a(j) .lt. a(j + 1) ) j = j + 1
1278 if (aa .lt. a(j))
then
1297 integer,
intent(in) :: n
1298 real(kind=
rp),
intent(inout) :: b(n)
1299 integer,
intent(in) :: ind(n)
1300 real(kind=
rp) :: temp(n)
1317 integer,
intent(in) :: n
1318 integer(i4),
intent(inout) :: b(n)
1319 integer,
intent(in) :: ind(n)
1320 integer(i4) :: temp(n)
1337 integer,
intent(in) :: n
1338 real(kind=
rp),
intent(inout) :: b(n)
1339 integer,
intent(in) :: ind(n)
1340 real(kind=
rp) :: temp(n)
1357 integer,
intent(in) :: n
1358 integer(i4),
intent(inout) :: b(n)
1359 integer,
intent(in) :: ind(n)
1360 integer(i4) :: temp(n)
1377 integer,
intent(in) :: n
1378 real(kind=
rp),
intent(inout) :: b(n)
1379 integer,
intent(inout) :: ind(n)
1380 real(kind=
rp) :: temp(n)
1381 integer :: tempind(n)
1387 tempind(jj) = ind(i)
1400 integer,
intent(in) :: n
1401 integer(i4),
intent(inout) :: b(n)
1402 integer,
intent(inout) :: ind(n)
1403 integer(i4) :: temp(n)
1404 integer :: tempind(n)
1410 tempind(jj) = ind(i)
1422 integer,
intent(in) :: n
1423 real(kind=
rp),
dimension(n),
intent(inout) :: a
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) =
max(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) =
max(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
1484 integer,
intent(in) :: n
1485 real(kind=
rp),
dimension(n),
intent(inout) :: a
1486 real(kind=
rp),
dimension(n),
intent(in) :: b
1490 a(i) = min(a(i), b(i))
1496 integer,
intent(in) :: n
1497 real(kind=
rp),
dimension(n),
intent(inout) :: a
1498 real(kind=
rp),
dimension(n),
intent(in) :: b, c
1502 a(i) = min(b(i), c(i))
1508 integer,
intent(in) :: n
1509 real(kind=
rp),
dimension(n),
intent(inout) :: a
1510 real(kind=
rp),
intent(in) :: b
1520 integer,
intent(in) :: n
1521 real(kind=
rp),
dimension(n),
intent(inout) :: a
1522 real(kind=
rp),
dimension(n),
intent(in) :: b
1523 real(kind=
rp),
intent(in) :: c
1534 function matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33) &
1536 real(kind=
rp),
intent(in) :: a11, a12, a13, a21, a22, a23, a31, a32, a33
1557 real(kind=
xp),
intent(in) :: a(3,3)
1558 real(kind=
xp) :: b(3,3)
1559 real(kind=
xp) :: detinv
1563 detinv = 1.0_xp /
real(a(1,1)*a(2,2)*a(3,3) - a(1,1)*a(2,3)*a(3,2) &
1564 - a(1,2)*a(2,1)*a(3,3) + a(1,2)*a(2,3)*a(3,1)&
1565 + a(1,3)*a(2,1)*a(3,2) - a(1,3)*a(2,2)*a(3,1),
xp)
1568 b(1,1) = +detinv * (a(2,2)*a(3,3) - a(2,3)*a(3,2))
1569 b(2,1) = -detinv * (a(2,1)*a(3,3) - a(2,3)*a(3,1))
1570 b(3,1) = +detinv * (a(2,1)*a(3,2) - a(2,2)*a(3,1))
1571 b(1,2) = -detinv * (a(1,2)*a(3,3) - a(1,3)*a(3,2))
1572 b(2,2) = +detinv * (a(1,1)*a(3,3) - a(1,3)*a(3,1))
1573 b(3,2) = -detinv * (a(1,1)*a(3,2) - a(1,2)*a(3,1))
1574 b(1,3) = +detinv * (a(1,2)*a(2,3) - a(1,3)*a(2,2))
1575 b(2,3) = -detinv * (a(1,1)*a(2,3) - a(1,3)*a(2,1))
1576 b(3,3) = +detinv * (a(1,1)*a(2,2) - a(1,2)*a(2,1))
1582 real(kind=
rp),
intent(in) :: x
1583 real(kind=
rp) :: val
1584 real(kind=
rp),
parameter :: xdmin = 0.0001_rp
1585 real(kind=
rp),
parameter :: xdmax = 0.9999_rp
1588 if (x <= xdmin)
then
1591 else if (x >= xdmax)
then
1596 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1599 val = 1.0_rp / (1.0_rp + exp(g))
1605 real(kind=
rp),
intent(in) :: x
1606 real(kind=
rp) :: val
1607 real(kind=
rp),
parameter :: xdmin = 0.0001_rp
1608 real(kind=
rp),
parameter :: xdmax = 0.9999_rp
1609 real(kind=
rp) :: arg, g, dg, s_val
1611 if (x <= xdmin .or. x >= xdmax)
then
1618 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1621 dg = -(1.0_rp / ((x - 1.0_rp)**2)) - (1.0_rp / (x**2))
1624 s_val = 1.0_rp / (1.0_rp + exp(g))
1626 val = -s_val * (1.0_rp - s_val) * dg
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
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 .
subroutine, public face_masked_gather_copy_0(a, b, mask, facet, lx, ly, lz, n_mask)
Gather values from a face-local SEM field to a reduced contiguous 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 .
subroutine, public cwrap(a, min_val, max_val, n)
Wrap value around a range [min, max)
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.