25 use,
intrinsic :: iso_c_binding, only: c_double
26 use,
intrinsic :: iso_fortran_env
112 real(c_double) function expm1(x) bind(c, name=
'expm1')
114 real(c_double),
intent(in),
value :: x
139 real(real64),
intent(in) :: x, y
140 real(real64),
optional,
intent(in) :: rtol
141 real(real64),
optional,
intent(in) :: atol
142 real(real64) :: atol_, rtol_
144 atol_ = optional_default(atol, 1.e-8_real64)
145 rtol_ = optional_default(rtol, 1.e-5_real64)
152 complex(real64),
intent(in) :: x, y
153 real(real64),
optional,
intent(in) :: rtol
154 real(real64),
optional,
intent(in) :: atol
155 real(real64) :: atol_, rtol_
157 atol_ = optional_default(atol, 1.e-8_real64)
158 rtol_ = optional_default(rtol, 1.e-5_real64)
168 integer,
intent(in) :: dim
169 integer,
intent(in) :: diag
170 integer :: matrix(dim, dim)
176 matrix(ii, ii) = diag
182 recursive function hermite(n, x)
result (h)
183 integer,
intent(in) :: n
184 real(real64),
intent(in) :: x
203 integer,
intent(in) :: n
219 real(real64),
intent(in) :: xx(3)
220 integer,
intent(in) :: li, mi
221 complex(real64),
intent(out) :: ylm
224 real(real64) :: dx, dy, dz, r, plm, cosm, sinm, cosmm1, sinmm1, cosphi, sinphi
225 real(real64),
parameter :: tiny = 1.e-15_real64
231 ylm = cmplx(
sqrt(m_one/(m_four*m_pi)), m_zero, real64)
235 r =
sqrt(xx(1)**2 + xx(2)**2 + xx(3)**2)
248 if(abs(dz) > m_one) dz = sign(m_one, dz)
249 plm = loct_legendre_sphplm(li, abs(mi), dz)
252 r =
sqrt(dx*dx + dy*dy)
266 cosm = cosmm1*cosphi - sinmm1*sinphi
267 sinm = cosmm1*sinphi + sinmm1*cosphi
271 ylm = plm*cmplx(cosm, sinm, real64)
290 real(real64),
intent(in) :: xx(3)
291 integer,
intent(in) :: li, mi
292 real(real64),
intent(out) :: ylm
294 integer,
parameter :: lmaxd = 20
295 real(real64),
parameter :: tiny = 1.e-15_real64
296 integer :: i, ilm0, l, m, mabs
297 integer,
save :: lmax = -1
299 real(real64) :: cmi, cosm, cosmm1, cosphi, dphase, dplg, fac, &
300 fourpi, plgndr, phase, pll, pmm, pmmp1, sinm, &
301 sinmm1, sinphi, rsize, rx, ry, rz, xysize
302 real(real64),
save :: c(0:(lmaxd+1)*(lmaxd+1))
313 do i = l - m + 1, l + m
316 c(ilm0 + m) =
sqrt(fac)
318 if (m /= 0) c(ilm0 + m) = c(ilm0 + m)*
sqrt(m_two)
319 c(ilm0 - m) = c(ilm0 + m)
332 rsize =
sqrt(xx(1)**2 + xx(2)**2 + xx(3)**2)
333 if (rsize < tiny)
then
358 ylm = c(4)*6.0_real64*rx*ry
360 ylm = (-c(5))*m_three*ry*rz
362 ylm = c(6)*m_half*(m_three*rz*rz - m_one)
364 ylm = (-c(7))*m_three*rx*rz
366 ylm = c(8)*m_three*(rx*rx - ry*ry)
373 xysize =
sqrt(max(rx*rx + ry*ry, tiny))
381 cosm = cosmm1*cosphi - sinmm1*sinphi
382 sinm = cosmm1*sinphi + sinmm1*cosphi
390 dphase = (-mabs)*sinm
396 if (mabs > m_zero)
then
398 pmm = (-pmm)*fac*xysize
405 dplg = (-li)*rz*pmm/(xysize**2)
407 pmmp1 = rz*(2*mabs + 1)*pmm
408 if (li == mabs + 1)
then
410 dplg = -((li*rz*pmmp1 - (mabs + li)*pmm)/(xysize**2))
413 pll = (rz*(2*l - 1)*pmmp1 - (l + mabs - 1)*pmm)/(l - mabs)
418 dplg = -((li*rz*pll - (l + mabs - 1)*pmm)/(xysize**2))
424 ylm = cmi*plgndr*phase
439 subroutine weights(N, M, cc, side)
440 integer,
intent(in) :: n, m
441 real(real64),
intent(out) :: cc(0:,0:,0:)
442 integer,
optional,
intent(in) :: side
444 integer :: i, j, k, mn, side_
445 real(real64) :: c1, c2, c3, c4, c5, xi
446 real(real64),
allocatable :: x(:)
450 safe_allocate(x(0:m))
452 if (
present(side))
then
462 x(:) = (/(-i,i=0,mn)/)
466 x(:) = (/(i,i=0,mn)/)
470 x(:) = (/0,(-i,i,i=1,mn)/)
493 if (j <= n) cc(k, j - 1, j) = m_zero
494 cc(k, j, 0) = c4*cc(k, j - 1, 0)/c3
497 cc(k, j, i) = (c4*cc(k, j - 1, i) - i*cc(k, j - 1, i - 1))/c3
500 cc(j, j, 0) = -c1*c5*cc(j - 1, j - 1, 0) / c2
504 cc(j, j, i) = c1*(i*cc(j - 1, j - 1, i - 1) - c5*cc(j - 1, j - 1, i))/c2
515 real(real64)
pure function
ddelta(i, j)
516 integer,
intent(in) :: i
517 integer,
intent(in) :: j
531 integer,
intent(in) :: nn
532 real(real64),
intent(in) :: xa(:)
533 real(real64),
intent(in) :: xx
534 real(real64),
intent(out) :: cc(:)
544 cc(ii) = cc(ii)*(xx - xa(kk))/(xa(ii) - xa(kk))
553 logical pure function
even(n)
554 integer,
intent(in) :: n
556 even = (mod(n, 2) == 0)
563 logical pure function
odd(n)
564 integer,
intent(in) :: n
574 real(real64),
intent(in) :: x(:)
575 real(real64),
intent(out) :: u(:)
578 real(real64) :: sumx2
579 real(real64),
allocatable :: xx(:)
585 assert(
size(u) == n-1)
588 safe_allocate(xx(1:n))
590 if (abs(x(j))<1.0e-8_real64)
then
601 sumx2 = sumx2 + xx(j)**2
603 if (abs(sumx2) <= m_epsilon .and. abs(xx(k)) <= m_epsilon)
exit
607 u(n-1) =
atan2(xx(n), xx(n-1))
621 real(real64),
intent(in) :: u(:)
622 real(real64),
intent(out) :: x(:)
630 assert(
size(u) == n-1)
637 x(2) =
sin(u(1))*
cos(u(2))
638 x(3) =
sin(u(1))*
sin(u(2))
641 x(2) =
sin(u(1))*
cos(u(2))
646 x(j) = x(j) *
sin(u(k))
652 x(n) = x(n) *
sin(u(k))
654 x(n) = x(n) *
sin(u(n-1))
667 real(real64),
intent(out) :: grad_matrix(:,:)
668 real(real64),
intent(in) :: r
669 real(real64),
intent(in) :: x(:)
679 grad_matrix(1,1) = -r*
sin(x(1))
681 grad_matrix(m,1) = m_zero
688 grad_matrix(m,l) = -r*grad_matrix(m,l)*product(
sin(x(1:l)))
690 grad_matrix(m,l) = grad_matrix(m,l)*r*
cos(x(m))*
cos(x(l))*product(
sin(x(1:m-1)))*product(
sin(x(m+1:l-1)))
692 grad_matrix(m,l) = m_zero
699 grad_matrix(m,n) = r*
cos(x(m))*grad_matrix(m,n)*product(
sin(x(1:m-1)))*product(
sin(x(m+1:n-1)))
706 integer(int64) pure function
pad88(size, blk)
707 integer(int64),
intent(in) :: size
708 integer(int64),
intent(in) :: blk
716 pad88 =
size + blk - mm
721 integer(int64) pure function
pad48(size, blk)
722 integer,
intent(in) :: size
723 integer(int64),
intent(in) :: blk
729 integer(int64) pure function
pad8(size, blk)
730 integer(int64),
intent(in) :: size
731 integer,
intent(in) :: blk
737 integer pure function
pad4(size, blk)
738 integer,
intent(in) :: size
739 integer,
intent(in) :: blk
741 pad4 = int(
pad88(int(
size, int64), int(blk, int64)), int32)
750 integer pure function
pad_pow2(size)
751 integer,
intent(in) :: size
759 if (
size == 0)
return
765 if (mm - 2*mm2 /= 0)
exit
779 real(real64)
pure function
dlog2(xx)
780 real(real64),
intent(in) :: xx
787 integer pure function
ilog2(xx)
788 integer,
intent(in) :: xx
795 integer(int64) pure function
llog2(xx)
796 integer(int64),
intent(in) :: xx
798 llog2 = nint(
log2(real(xx, real64) ), kind=int64)
804 complex(real64),
intent(in) :: z
815 complex(real64) pure function
phi1(z)
816 complex(real64),
intent(in) :: z
818 real(real64),
parameter :: cut = 0.002_real64
820 if (abs(z) < cut)
then
821 phi1 = m_one + m_half*z*(m_one + z/m_three*(m_one + m_fourth*z*(m_one + z/m_five)))
834 complex(real64) pure function
phi2(z)
835 complex(real64),
intent(in) :: z
837 real(real64),
parameter :: cut = 0.002_real64
839 if (abs(z) < cut)
then
840 phi2 = m_half + z/6.0_real64*(m_one + m_fourth*z*(m_one + z/m_five*(m_one + z/6.0_real64)))
842 phi2 = (
exp(z) - z - m_z1)/z**2
849 real(real64),
intent(in) :: x
857 integer,
intent(in) :: n
864 message(1) =
"Internal error: is_prime does not take negative numbers."
865 call messages_fatal(1)
873 root = nint(
sqrt(real(n, real64) ))
875 if (mod(n,i) == 0)
then
892 real(real64),
intent(out) :: r(:,:)
893 real(real64),
intent(in) :: ff(:)
894 real(real64),
intent(in) :: tt(:)
897 real(real64) :: th, uv, uu, vv, ft
898 real(real64),
allocatable :: axis(:), u(:), v(:),
f(:), t(:), p(:)
904 assert((dim < 3) .or. (dim > 2))
905 assert(
size(tt,1) == dim)
906 assert((
size(r,1) == dim) .and. (
size(r,2) == dim))
908 safe_allocate(u(1:dim))
909 safe_allocate(v(1:dim))
910 safe_allocate(
f(1:dim))
911 safe_allocate(t(1:dim))
918 ft = dot_product(
f,t)
920 if (abs(ft) < m_one)
then
933 safe_allocate(axis(1:dim))
936 u =
f / dot_product(
f,
f)
937 v = t /dot_product(t,t)
940 axis = axis / norm2(axis)
942 r(1,1) =
cos(th) + axis(1)**2 * (1 -
cos(th))
943 r(1,2) = axis(1)*axis(2)*(1-
cos(th)) + axis(3)*
sin(th)
944 r(1,3) = axis(1)*axis(3)*(1-
cos(th)) - axis(2)*
sin(th)
946 r(2,1) = axis(2)*axis(1)*(1-
cos(th)) - axis(3)*
sin(th)
947 r(2,2) =
cos(th) + axis(2)**2 * (1 -
cos(th))
948 r(2,3) = axis(2)*axis(3)*(1-
cos(th)) + axis(1)*
sin(th)
950 r(3,1) = axis(3)*axis(1)*(1-
cos(th)) + axis(2)*
sin(th)
951 r(3,2) = axis(3)*axis(2)*(1-
cos(th)) - axis(1)*
sin(th)
952 r(3,3) =
cos(th) + axis(3)**2 * (1 -
cos(th))
954 safe_deallocate_a(axis)
963 r(1,1) = u(1)**2 + (1-u(1)**2)*
cos(th)
964 r(1,2) = u(1)*u(2)*(1-
cos(th)) - u(3)*
sin(th)
965 r(1,3) = u(1)*u(3) + u(2)*
sin(th)
967 r(2,1) = u(1)*u(2)*(1-
cos(th)) + u(3)*
sin(th)
968 r(2,2) = u(2)**2 + (1-u(2)**2)*
cos(th)
969 r(2,3) = u(2)*u(3)*(1-
cos(th)) - u(1)*
sin(th)
971 r(3,1) = u(1)*u(3)*(1-
cos(th)) - u(2)*
sin(th)
972 r(3,2) = u(2)*u(3)*(1-
cos(th)) + u(1)*
sin(th)
973 r(3,3) = u(3)**2 + (1-u(3)**2)*
cos(th)
978 safe_allocate(p(1:dim))
980 if (abs(
f(1)) <= abs(
f(2)) .and. abs(
f(1)) < abs(
f(3)))
then
981 p = (/m_one, m_zero, m_zero/)
982 else if (abs(
f(2)) < abs(
f(1)) .and. abs(
f(2)) <= abs(
f(3)))
then
983 p = (/m_zero, m_one, m_zero/)
984 else if (abs(
f(3)) <= abs(
f(1)) .and. abs(
f(3)) < abs(
f(2)))
then
985 p = (/m_zero, m_zero, m_one/)
991 uu = dot_product(u,u)
992 vv = dot_product(v,v)
993 uv = dot_product(u,v)
998 r(i,j) =
ddelta(i,j) - m_two * u(i)*u(j)/uu - m_two * v(i)*v(j)/vv &
999 + 4_real64*uv * v(i)*u(j) /(uu*vv)
1004 safe_deallocate_a(p)
1020 safe_deallocate_a(u)
1021 safe_deallocate_a(v)
1037 real(real64),
intent(in) :: x, h
1038 real(real64),
intent(out) :: res, err
1041 use,
intrinsic :: iso_fortran_env
1043 real(real64),
intent(in) :: x
1044 real(real64),
intent(inout) :: fx
1048 real(real64),
parameter :: con = 1.4_real64, big = 1.0e30_real64, safe = m_two
1049 integer,
parameter :: ntab = 20
1051 real(real64) :: errt, fac, hh, fx1, fx2
1052 real(real64),
allocatable :: a(:, :)
1056 if (abs(h) <= m_epsilon)
then
1057 message(1) =
"h must be nonzero in numder_ridders"
1058 call messages_fatal(1)
1061 safe_allocate(a(1:ntab, 1:ntab))
1066 a(1,1) = (fx1-fx2) / (m_two*hh)
1072 a(1,i) = (fx1-fx2) / (m_two*hh)
1075 a(j,i) = (a(j-1,i)*fac-a(j-1,i-1)) / (fac-m_one)
1077 errt = max(abs(a(j,i)-a(j-1,i)),abs(a(j,i)-a(j-1,i-1)))
1078 if (errt .le. err)
then
1083 if (abs(a(i,i)-a(i-1,i-1)) .ge. safe*err)
return
1086 safe_deallocate_a(a)
1092 real(real64),
intent(in) :: a(:)
1093 complex(real64),
intent(in) :: b(:)
1095 complex(real64) :: c(1:3)
1097 c(1) = a(2)*b(3) - a(3)*b(2)
1098 c(2) = a(3)*b(1) - a(1)*b(3)
1099 c(3) = a(1)*b(2) - a(2)*b(1)
1105 complex(real64),
intent(in) :: a(:)
1106 real(real64),
intent(in) :: b(:)
1108 complex(real64) :: c(1:3)
1110 c(1) = a(2)*b(3) - a(3)*b(2)
1111 c(2) = a(3)*b(1) - a(1)*b(3)
1112 c(3) = a(1)*b(2) - a(2)*b(1)
1218 integer,
intent(in) :: np
1219 integer,
intent(in) :: nn, mm
1220 real(real64),
intent(in) :: xx(np)
1221 real(real64),
intent(out) :: cx(np)
1224 real(real64),
allocatable :: cx_tmp(:,:)
1228 safe_allocate(cx_tmp(1:np, 0:nn))
1240 cx_tmp(1:np,0) = m_one
1241 cx_tmp(1:np,1) = real(mm + 1, real64) - xx(1:np)
1244 cx_tmp(1:np, ii) = &
1245 (( real(mm + 2*ii - 1, real64) - xx(1:np)) * cx_tmp(1:np,ii-1) &
1246 + real(-mm -ii + 1, real64) * cx_tmp(1:np,ii-2)) / real(ii, real64)
1249 cx(1:np) = cx_tmp(1:np, nn)
1251 safe_deallocate_a(cx_tmp)
1256 integer,
intent(in) :: nn
1257 real(real64),
intent(inout) :: aa(:, :)
1263 aa(jj, ii) = aa(ii, jj)
1269 integer,
intent(in) :: nn
1270 complex(real64),
intent(inout) :: aa(:, :)
1276 aa(jj, ii) = conjg(aa(ii, jj))
1278 aa(jj, jj) = real(aa(jj, jj), real64)
1283 integer,
intent(in) :: nn
1284 real(real64),
intent(inout) :: aa(:, :)
1290 aa(ii, jj) = aa(jj, ii)
1296 integer,
intent(in) :: nn
1297 complex(real64),
intent(inout) :: aa(:, :)
1303 aa(ii, jj) = conjg(aa(jj, ii))
1305 aa(jj, jj) = real(aa(jj, jj), real64)
1311 integer,
intent(in) :: nn
1312 real(real64),
intent(inout) :: aa(:, :)
1318 aa(ii, jj) = m_half * (aa(ii, jj) + aa(jj, ii))
1319 aa(jj, ii) = aa(ii, jj)
1326 integer,
intent(in) :: n
1327 complex(real64),
intent(inout) :: a(:, :)
1333 a(i, j) = m_half * (a(i, j) + conjg(a(j, i)))
1334 a(j, i) = conjg(a(i, j))
1341 real(real64),
intent(inout) :: aa(:, :)
1351 real(real64),
intent(in) :: a(:, :)
1352 real(real64),
optional,
intent(in) :: tol
1354 real(real64) :: tolerance
1358 tolerance = optional_default(tol, 1.0e-8_real64)
1360 assert(
size(a, 1) ==
size(a, 2))
1369 logical pure function is_diagonal(dim, matrix, tol)
1370 integer,
intent(in) :: dim
1371 real(real64),
intent(in) :: matrix(:, :)
1372 real(real64),
intent(in) :: tol
1375 real(real64) :: max_diag
1379 max_diag = max(abs(matrix(ii, ii)), max_diag)
1386 if (abs(matrix(ii, jj)) > tol*max_diag)
is_diagonal = .false.
1394 integer,
intent(in) :: num
1395 integer,
intent(in) :: base
1396 integer,
intent(out) :: converted(:)
1398 integer :: i, num_, length
1406 assert(
size(converted) >= length)
1410 do i = 1,
size(converted)
1411 converted(i) = mod(num_, base)
1420 integer,
intent(in) :: converted(:)
1421 integer,
intent(in) :: base
1422 integer,
intent(out) :: num
1429 do i = 1,
size(converted)
1430 num = num + converted(i) * base**(i-1)
1436#include "complex.F90"
1437#include "math_inc.F90"
1441#include "math_inc.F90"
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
double acos(double __x) __attribute__((__nothrow__
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
double atan2(double __y, double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
This module is intended to contain "only mathematical" functions and procedures.
subroutine zsymmetrize_matrix(n, a)
make a complex matrix Hermitian
complex(real64) pure function, public phi2(z)
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2.
integer pure function ilog2(xx)
logical pure function, public even(n)
Returns if n is even.
pure real(real64) function, dimension(dim, dim) ddiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer versi...
subroutine, public generalized_laguerre_polynomial(np, nn, mm, xx, cx)
complex(real64) pure function, public exponential(z)
Wrapper for exponential.
subroutine dinterpolate_0(xa, ya, x, y)
subroutine, public ylmr_cmplx(xx, li, mi, ylm)
Computes spherical harmonics ylm at position (x, y, z)
subroutine zinterpolate_0(xa, ya, x, y)
subroutine, public convert_to_base(num, base, converted)
convert an integer in base 10 to a different base
subroutine dinterpolate_1(xa, ya, x, y)
real(real64) pure function dlog2(xx)
integer(int64) pure function pad88(size, blk)
real(real64) pure function, public square_root(x)
Wrapper for sqrt.
integer pure function pad4(size, blk)
subroutine dlower_triangular_to_hermitian(nn, aa)
pure complex(real64) function, dimension(1:3), public dzcross_product(a, b)
subroutine, public cartesian2hyperspherical(x, u)
Performs a transformation of an n-dimensional vector from Cartesian coordinates to hyperspherical coo...
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
subroutine, public dzero_small_elements_matrix(aa, tol)
subroutine zinterpolate_1(xa, ya, x, y)
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
integer(int64) pure function pad48(size, blk)
subroutine dsymmetrize_matrix(nn, aa)
symmetrize a real matrix
subroutine, public ylmr_real(xx, li, mi, ylm)
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x,...
subroutine dupper_triangular_to_hermitian(nn, aa)
subroutine, public interpolation_coefficients(nn, xa, xx, cc)
logical function, public is_symmetric(a, tol)
Check if a 2D array is symmetric.
subroutine, public generate_rotation_matrix(R, ff, tt)
Generates a rotation matrix R to rotate a vector f to t.
elemental logical function zis_close_scalar(x, y, rtol, atol)
Same as dis_close_scalar for complex numbers.
logical pure function, public odd(n)
Returns if n is odd.
real(real64) pure function, public ddelta(i, j)
logical function, public is_prime(n)
complex(real64) pure function, public phi1(z)
Compute phi1(z) = (exp(z)-1)/z.
pure complex(real64) function, dimension(1:3), public zcross_product(a, b)
subroutine, public weights(N, M, cc, side)
Compute the weights for finite-difference calculations:
subroutine zinterpolate_2(xa, ya, x, y)
pure complex(real64) function, dimension(dim, dim) zdiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer versi...
subroutine zlower_triangular_to_hermitian(nn, aa)
logical pure function, public is_diagonal(dim, matrix, tol)
Returns true is the matrix of size dim x dim is diagonal, given a relative tolerance.
subroutine zupper_triangular_to_hermitian(nn, aa)
pure integer function, dimension(dim, dim) idiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the real and comp...
integer(int64) pure function pad8(size, blk)
subroutine dinterpolate_2(xa, ya, x, y)
subroutine, public numder_ridders(x, h, res, err, f)
Numerical derivative (Ridder`s algorithm).
integer(int64) pure function llog2(xx)
recursive real(real64) function, public hermite(n, x)
elemental logical function dis_close_scalar(x, y, rtol, atol)
Are and equal within a tolerance.
pure complex(real64) function, dimension(1:3), public zdcross_product(a, b)
subroutine, public convert_from_base(converted, base, num)
convert an integer to base 10 from a different base
subroutine, public hyperspherical2cartesian(u, x)
Performs the inverse transformation of cartesian2hyperspherical.
subroutine, public hypersphere_grad_matrix(grad_matrix, r, x)
Gives the hyperspherical gradient matrix, which contains the derivatives of the Cartesian coordinates...
recursive integer function, public factorial(n)
static double f(double w, void *p)