24 use,
intrinsic :: iso_fortran_env
31#ifdef HAVE_LIBXC_FUNCS
51 integer,
parameter :: VALCONF_STRING_LENGTH = 80
56 character(len=3) :: symbol =
""
61 real(real64) :: occ(15,2) =
m_zero
62 real(real64) :: j(15) =
m_zero
70 type(valconf_t),
intent(out) :: cout
71 type(valconf_t),
intent(in) :: cin
76 cout%symbol = cin%symbol
88 type(valconf_t),
intent(inout) :: c
89 integer,
intent(in) :: lmax
92 integer,
allocatable :: order(:)
93 type(valconf_t) :: tmp_c
100 safe_allocate(order(1:lmax+1))
101 call sort(c%l(1:lmax+1), order)
103 c%occ(ll+1, 1:2) = tmp_c%occ(order(ll+1), 1:2)
104 c%n(ll+1) = tmp_c%n(order(ll+1))
105 c%j(ll+1) = tmp_c%j(order(ll+1))
107 safe_deallocate_a(order)
114 type(valconf_t),
intent(inout) :: conf
124 conf%occ(ii, 1) = min(xx, real(2*ll+1, real64) )
125 conf%occ(ii, 2) = xx - conf%occ(ii,1)
133 subroutine atomhxc(namespace, functl, g, nspin, dens, v, extra)
134 type(namespace_t),
intent(in) :: namespace
135 character(len=*),
intent(in) :: functl
136 type(logrid_t),
intent(in) :: g
137 integer,
intent(in) :: nspin
138 real(real64),
intent(in) :: dens(g%nrval, nspin)
139 real(real64),
intent(out) :: v(g%nrval, nspin)
140 real(real64),
intent(in),
optional :: extra(g%nrval)
142 character(len=5) :: xcfunc, xcauth
144 real(real64),
allocatable :: xc(:, :), ve(:, :), rho(:, :)
145 real(real64) :: r2, ex, ec, dx, dc
149 safe_allocate( ve(1:g%nrval, 1:nspin))
150 safe_allocate( xc(1:g%nrval, 1:nspin))
151 safe_allocate(rho(1:g%nrval, 1:nspin))
158 rho(:, 1) = rho(:, 1) + dens(:, is)
162 call vhrtre(rho(:, 1), ve(:, is), g%rofi, g%drdi, g%s, g%nrval, g%a)
164 ve(1, 1:nspin) = ve(2, 1:nspin)
174 message(1) =
'Internal Error in atomhxc: unknown functl'
178 rho(:, :) = dens(:, :)
180 if (
present(extra)) rho(:, is) = rho(:, is) + extra(:)/nspin
181 rho(2:g%nrval, is) = rho(2:g%nrval, is)/(
m_four*
m_pi*g%r2ofi(2:g%nrval))
183 r2 = g%rofi(2)/(g%rofi(3)-g%rofi(2))
184 rho(1, 1:nspin) = rho(2, 1:nspin) - (rho(3, 1:nspin)-rho(2, 1:nspin))*r2
192 call atomxc(namespace, xcfunc, xcauth, g%nrval, g%nrval, g%rofi, nspin, rho, ex, ec, dx, dc, xc)
195 safe_deallocate_a(ve)
196 safe_deallocate_a(xc)
197 safe_deallocate_a(rho)
242 subroutine atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
244 character(len=*),
intent(in) :: functl, author
245 integer,
intent(in) :: nr, maxr
246 real(real64),
intent(in) :: rmesh(maxr)
247 integer,
intent(in) :: nspin
248 real(real64),
intent(in) :: dens(maxr,nspin)
249 real(real64),
intent(out) :: dc, dx, ec, ex, v_xc(maxr,nspin)
257 integer,
parameter :: nn = 5
265 real(real64),
parameter :: eunit =
m_half
271 real(real64),
parameter :: dvmin = 1.0e-12_real64
278 integer :: in, in1, in2, ir, is, jn
279 real(real64) :: aux(maxr), d(nspin), &
280 dgdm(-nn:nn), dgidfj(-nn:nn), drdm, dvol, &
281 f1, f2, gd(3,nspin), sigma(3), vxsigma(3), vcsigma(3), &
283 dexdd(nspin), dexdgd(3,nspin), decdd(nspin), decdgd(3,nspin)
284 type(xc_f03_func_t) :: x_conf, c_conf
289 assert(nspin == 1 .or. nspin == 2)
292 if (functl .eq.
'LDA' .or. functl .eq.
'lda' .or. &
293 functl .eq.
'LSD' .or. functl .eq.
'lsd')
then
295 else if (functl .eq.
'GGA' .or. functl .eq.
'gga')
then
299 write(
message(1),
'(a,a)')
'atomxc: Unknown functional ', functl
305 call xc_f03_func_init(x_conf, xc_gga_x_pbe, nspin)
306 call xc_f03_func_init(c_conf, xc_gga_c_pbe, nspin)
308 call xc_f03_func_init(x_conf, xc_lda_x, nspin)
309 if (author .eq.
'CA' .or. author .eq.
'ca' .or. &
310 author .eq.
'PZ' .or. author .eq.
'pz')
then
311 call xc_f03_func_init(c_conf, xc_lda_c_pz, nspin)
312 else if (author .eq.
'PW92' .or. author .eq.
'pw92')
then
313 call xc_f03_func_init(c_conf, xc_lda_c_pw, nspin)
315 write(
message(1),
'(a,a)')
'LDAXC: Unknown author ', author
341 in1 = max( 1, ir-nn ) - ir
342 in2 = min( nr, ir+nn ) - ir
350 if (jn /= 0) dgdm(in) = dgdm(in) +
m_one / (0 - jn)
356 if (jn /= in .and. jn /= 0) f1 = f1 * (0 - jn)
357 if (jn /= in) f2 = f2 * (in - jn)
366 drdm = drdm + rmesh(ir+in) * dgdm(in)
370 dvol = 4 *
m_pi * rmesh(ir)**2 * drdm
373 dvol = safe_tol(dvol, dvmin)
374 if (ir == 1 .or. ir == nr) dvol = dvol / 2
375 if (gga) aux(ir) = dvol
381 dgidfj(in) = dgdm(in) / drdm
393 gd(3,is) = gd(3,is) + dgidfj(in) * dens(ir+in,is)
402 sigma(1) = gd(3, 1)*gd(3, 1)
404 sigma(2) = gd(3, 1) * gd(3, 2)
405 sigma(3) = gd(3, 2) * gd(3, 2)
411 call xc_f03_gga_exc_vxc(x_conf, int(1, c_size_t), d, sigma, epsx, dexdd, vxsigma)
412 call xc_f03_gga_exc_vxc(c_conf, int(1, c_size_t), d, sigma, epsc, decdd, vcsigma)
414 call xc_f03_lda_exc_vxc(x_conf, int(1, c_size_t), d, epsx, dexdd)
415 call xc_f03_lda_exc_vxc(c_conf, int(1, c_size_t), d, epsc, decdd)
422 dexdgd(:, 1) =
m_two * vxsigma(1) * gd(:, 1)
423 decdgd(:, 1) =
m_two * vcsigma(1) * gd(:, 1)
425 dexdgd(:, 1) = dexdgd(:, 1) + vxsigma(2)*gd(:, 2)
426 dexdgd(:, 2) =
m_two * vxsigma(3) * gd(:, 2) + vxsigma(2) * gd(:, 1)
427 decdgd(:, 1) = decdgd(:, 1) + vcsigma(2)*gd(:, 2)
428 decdgd(:, 2) =
m_two * vcsigma(3) * gd(:, 2) + vcsigma(2) * gd(:, 1)
435 ex = ex + dvol * d(is) * epsx(1)
436 ec = ec + dvol * d(is) * epsc(1)
437 dx = dx + dvol * d(is) * (epsx(1) - dexdd(is))
438 dc = dc + dvol * d(is) * (epsc(1) - decdd(is))
440 v_xc(ir,is) = v_xc(ir,is) + dvol * (dexdd(is) + decdd(is))
442 dx = dx - dvol * dens(ir+in,is) * dexdgd(3,is) * dgidfj(in)
443 dc = dc - dvol * dens(ir+in,is) * decdgd(3,is) * dgidfj(in)
444 v_xc(ir+in,is) = v_xc(ir+in,is) + dvol * &
445 (dexdgd(3,is) + decdgd(3,is)) * dgidfj(in)
448 v_xc(ir,is) = dexdd(is) + decdd(is)
462 v_xc(ir,is) = v_xc(ir,is) / dvol
477 v_xc(ir,is) = v_xc(ir,is) / eunit
481 call xc_f03_func_end(x_conf)
482 call xc_f03_func_end(c_conf)
507 subroutine vhrtre(rho, v, r, drdi, srdrdi, nr, a)
508 real(real64),
intent(in) :: rho(:)
509 real(real64),
intent(out) :: v(:)
510 real(real64),
intent(in) :: r(:), drdi(:), srdrdi(:)
511 integer,
intent(in) :: nr
512 real(real64),
intent(in) :: a
514 real(real64) :: x,y,q,a2by4,ybyq,qbyy,qpartc,v0,qt,dz,t,beta,dv
515 integer :: nrm1,nrm2,ir
522 ybyq =
m_one - a*a/48.0_real64
537 dz = drdi(ir)*rho(ir)
543 dz = drdi(ir)*rho(ir)
547 dz = drdi(nr)*rho(nr)
549 qt = (qt + qt + dz)*
m_half
550 v0 = (v0 + v0 + dz/r(nr))
561 beta = drdi(ir)*t*rho(ir)
579 q = (y - beta / 12.0_real64) * qbyy
580 v(ir) =
m_two * t * q
587 x = x + a2by4 * q - beta
589 t = srdrdi(ir) / r(ir)
590 beta = t * drdi(ir) * rho(ir)
591 q = (y - beta / 12.0_real64) * qbyy
592 v(ir) =
m_two * t * q
604 qpartc = r(nr)*v(nr)/
m_two
645 subroutine egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
647 real(real64),
intent(in) :: h(:), s(:)
648 integer,
intent(in) :: n
649 real(real64),
intent(inout) :: e
650 real(real64),
intent(out) :: g(:), y(:)
651 integer,
intent(in) :: l
652 real(real64),
intent(in) :: z, a, b, rmax
653 integer,
intent(in) :: nprin, nnode
654 real(real64),
intent(in) :: dr
655 integer,
intent(out) :: ierr
657 integer :: i,ncor,n1,n2,niter,nt
658 real(real64) :: e1, e2, del, de, et, t
659 real(real64),
parameter :: tol = 1.0e-5_real64
707 if (.not. (et <= e1 .or. et >= e2 .or. nt < nnode-1 .or. nt > nnode))
then
709 if (abs(de) < tol)
then
714 call yofe(e,de,dr,rmax,h,s,y,n,l,ncor,nt,z,a,b)
737 if (n2 >= nnode) cycle
757 if (n1 < nnode) cycle
776 g(i) = y(i) / (
m_one - t / 12._real64)
778 call nrmlzg(namespace, g, s, n)
785 subroutine yofe(e,de,dr,rmax,h,s,y,nmax,l,ncor,nnode,z,a,b)
806 real(real64),
intent(inout) :: e
807 real(real64),
intent(out) :: de
808 real(real64),
intent(in) :: dr, rmax
809 real(real64),
intent(in) :: h(:), s(:)
810 real(real64),
intent(out) :: y(:)
811 integer,
intent(in) :: nmax, l
812 integer,
intent(in) :: ncor
813 integer,
intent(out) :: nnode
814 real(real64),
intent(in) :: z, a, b
816 real(real64) :: y2, g, gsg, x, gin, gsgin, xin
817 integer :: n,knk,nndin,i
818 real(real64) :: zdr, yn, ratio, t
825 if (h(n)-e*s(n) <
m_one)
then
832 call bcorgn(e,h,s,l,zdr,y2)
841 call numout(e,h,s,y,ncor,knk,nnode,y2,g,gsg,x)
855 if (.not. (n < nmax .or. abs(dr) > 1.e3_real64))
then
856 call bcrmax(e,dr,rmax,h,s,n,yn,a,b)
870 call numin(e,h,s,y,n,nndin,yn,gin,gsgin,xin,knk)
884 gsg = gsg + gsgin * ratio * ratio
885 t = h(knk) - e * s(knk)
886 de = g * (x + xin + t * g) / gsg
887 nnode = nnode + nndin
888 if (de <
m_zero) nnode = nnode + 1
897 subroutine nrmlzg(namespace, g, s, n)
920 real(real64),
intent(inout) :: g(:)
921 real(real64),
intent(in) :: s(:)
922 integer,
intent(in) :: n
925 real(real64) :: norm, srnrm
934 if (mod(n,2) /= 1)
then
935 write(
message(1),
'(a,i6)')
' nrmlzg: n should be odd. n =', n
942 norm=norm + g(i)*s(i)*g(i)
947 norm=norm + g(i)*s(i)*g(i)
952 norm = norm + g(i) * s(i) * g(i)
964 subroutine bcorgn(e,h,s,l,zdr,y2)
965 real(real64),
intent(in) :: e
966 real(real64),
intent(in) :: h(:), s(:)
967 integer,
intent(in) :: l
968 real(real64),
intent(in) :: zdr
969 real(real64),
intent(out) :: y2
971 real(real64) :: t2, t3, d2, c0, c1, c2
981 d2 = -((24._real64 + 10.0_real64 * t2) / (12._real64 - t2))
997 c0=c0/(
m_one-0.75_real64*zdr)
1003 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1005 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1006 d2=(d2-c1)/(
m_one-c2)
1013 c0=c0/(
m_one-0.75_real64*zdr)
1018 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1020 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1021 d2=(d2-c1)/(
m_one-c2)
1031 subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
1032 real(real64),
intent(in) :: e, dr, rmax
1033 real(real64),
intent(in) :: h(:), s(:)
1034 integer,
intent(in) :: n
1035 real(real64),
intent(out) :: yn
1036 real(real64),
intent(in) :: a, b
1038 real(real64) :: tnm1, tn, tnp1, beta, dg, c1, c2, c3, dn
1042 tnm1=h(n-1)-e*s(n-1)
1044 tnp1=h(n+1)-e*s(n+1)
1048 c2=24._real64*dg/(12._real64-tn)
1049 dn=-((24._real64+10._real64*tn)/(12._real64-tn))
1051 c1= (
m_one-tnm1/6.0_real64)/(
m_one-tnm1/12._real64)
1052 c3=-((
m_one-tnp1/6.0_real64)/(
m_one-tnp1/12._real64))
1053 yn=-((
m_one-c1/c3)/(dn-c2/c3))
1059 subroutine numin(e,h,s,y,n,nnode,yn,g,gsg,x,knk)
1060 real(real64),
intent(in) :: e, h(:), s(:)
1061 real(real64),
intent(inout) :: y(:)
1062 integer,
intent(in) :: n
1063 integer,
intent(out) :: nnode
1064 real(real64),
intent(in) :: yn
1065 real(real64),
intent(out) :: g, gsg, x
1066 integer,
intent(in) :: knk
1075 g=y(n)/(
m_one-t/12._real64)
1080 g=y(i)/(
m_one-t/12._real64)
1094 do i = n - 2, knk, -1
1101 g = y(i) / (
m_one - t / 12._real64)
1102 gsg = gsg + g * s(i) * g
1116 end subroutine numin
1119 subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
1120 real(real64),
intent(in) :: e, h(:), s(:)
1121 real(real64),
intent(inout) :: y(:)
1122 integer,
intent(in) :: ncor
1123 integer,
intent(inout) :: knk
1124 integer,
intent(out) :: nnode
1125 real(real64),
intent(in) :: y2
1126 real(real64),
intent(out) :: g, gsg, x
1129 real(real64) :: t, xl
1136 g=y(2)/(
m_one-t/12._real64)
1140 g=y(3)/(
m_one-t/12._real64)
1161 g = y(i) / (
m_one - t / 12._real64)
1162 gsg = gsg + g * s(i) * g
1163 if (.not. (nnode < ncor .or. xl*x >
m_zero))
then
This is the common interface to a sorting routine. It performs the shell algorithm,...
double sqrt(double __x) __attribute__((__nothrow__
subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
subroutine, public atomhxc(namespace, functl, g, nspin, dens, v, extra)
subroutine nrmlzg(namespace, g, s, n)
subroutine yofe(e, de, dr, rmax, h, s, y, nmax, l, ncor, nnode, z, a, b)
subroutine, public atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
Finds total exchange-correlation energy and potential for a ! spherical electron density distribution...
subroutine, public valconf_unpolarized_to_polarized(conf)
integer, parameter, public valconf_string_length
Following stuff is for the "val_conf_t" data type, made to handle atomic configurations.
subroutine, public valconf_copy(cout, cin)
subroutine, public valconf_sort_by_l(c, lmax)
subroutine numin(e, h, s, y, n, nnode, yn, g, gsg, x, knk)
subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
subroutine bcorgn(e, h, s, l, zdr, y2)
subroutine, public vhrtre(rho, v, r, drdi, srdrdi, nr, a)
VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY....
subroutine, public egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
egofv determines the eigenenergy and wavefunction corresponding ! to a particular l,...
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
This module is intended to contain "only mathematical" functions and procedures.