24 use,
intrinsic :: iso_fortran_env
30#ifdef HAVE_LIBXC_FUNCS
51 character(len=1),
parameter :: &
52 spec_notation(0:3) = (/
's',
'p',
'd',
'f' /)
54 integer,
parameter :: VALCONF_STRING_LENGTH = 80
59 character(len=3) :: symbol =
""
64 real(real64) :: occ(15,2) =
m_zero
65 real(real64) :: j(15) =
m_zero
73 type(valconf_t),
intent(out) :: cout
74 type(valconf_t),
intent(in) :: cin
79 cout%symbol = cin%symbol
93 type(valconf_t),
intent(in) :: c
94 character(len=VALCONF_STRING_LENGTH),
intent(out) :: s
100 write(s,
'(i2,1x,a2,i1,1x,i1,a1,6(i1,a1,f6.3,a1))') c%z, c%symbol, c%type, c%p,
':',&
101 (c%n(j),spec_notation(c%l(j)),c%occ(j,1),
',',j=1,c%p)
109 type(namespace_t),
intent(in) :: namespace
110 character(len=VALCONF_STRING_LENGTH),
intent(in) :: s
111 type(valconf_t),
intent(out) :: c
114 character(len=1) :: lvalues(1:6)
118 read (s,
'(i2,1x,a2,i1,1x,i1,1x,6(i1,a1,f6.3,1x))') c%z, c%symbol, c%type, c%p,&
119 (c%n(j),lvalues(j),c%occ(j,1),j=1,c%p)
121 select case (lvalues(j))
143 type(valconf_t),
intent(inout) :: conf
153 conf%occ(ii, 1) = min(xx, real(2*ll+1, real64) )
154 conf%occ(ii, 2) = xx - conf%occ(ii,1)
162 subroutine atomhxc(namespace, functl, g, nspin, dens, v, extra)
164 character(len=*),
intent(in) :: functl
166 integer,
intent(in) :: nspin
167 real(real64),
intent(in) :: dens(g%nrval, nspin)
168 real(real64),
intent(out) :: v(g%nrval, nspin)
169 real(real64),
intent(in),
optional :: extra(g%nrval)
171 character(len=5) :: xcfunc, xcauth
173 real(real64),
allocatable :: xc(:, :), ve(:, :), rho(:, :)
174 real(real64) :: r2, ex, ec, dx, dc
178 safe_allocate( ve(1:g%nrval, 1:nspin))
179 safe_allocate( xc(1:g%nrval, 1:nspin))
180 safe_allocate(rho(1:g%nrval, 1:nspin))
187 rho(:, 1) = rho(:, 1) + dens(:, is)
191 call vhrtre(rho(:, 1), ve(:, is), g%rofi, g%drdi, g%s, g%nrval, g%a)
193 ve(1, 1:nspin) = ve(2, 1:nspin)
203 message(1) =
'Internal Error in atomhxc: unknown functl'
207 rho(:, :) = dens(:, :)
209 if (
present(extra)) rho(:, is) = rho(:, is) + extra(:)/nspin
210 rho(2:g%nrval, is) = rho(2:g%nrval, is)/(
m_four*
m_pi*g%r2ofi(2:g%nrval))
212 r2 = g%rofi(2)/(g%rofi(3)-g%rofi(2))
213 rho(1, 1:nspin) = rho(2, 1:nspin) - (rho(3, 1:nspin)-rho(2, 1:nspin))*r2
221 call atomxc(namespace, xcfunc, xcauth, g%nrval, g%nrval, g%rofi, nspin, rho, ex, ec, dx, dc, xc)
224 safe_deallocate_a(ve)
225 safe_deallocate_a(xc)
226 safe_deallocate_a(rho)
271 subroutine atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
273 character(len=*),
intent(in) :: functl, author
274 integer,
intent(in) :: nr, maxr
275 real(real64),
intent(in) :: rmesh(maxr)
276 integer,
intent(in) :: nspin
277 real(real64),
intent(in) :: dens(maxr,nspin)
278 real(real64),
intent(out) :: dc, dx, ec, ex, v_xc(maxr,nspin)
286 integer,
parameter :: nn = 5
294 real(real64),
parameter :: eunit =
m_half
300 real(real64),
parameter :: dvmin = 1.0e-12_real64
307 integer :: in, in1, in2, ir, is, jn
308 real(real64) :: aux(maxr), d(nspin), &
309 dgdm(-nn:nn), dgidfj(-nn:nn), drdm, dvol, &
310 f1, f2, gd(3,nspin), sigma(3), vxsigma(3), vcsigma(3), &
312 dexdd(nspin), dexdgd(3,nspin), decdd(nspin), decdgd(3,nspin)
313 type(xc_f03_func_t) :: x_conf, c_conf
318 assert(nspin == 1 .or. nspin == 2)
321 if (functl .eq.
'LDA' .or. functl .eq.
'lda' .or. &
322 functl .eq.
'LSD' .or. functl .eq.
'lsd')
then
324 else if (functl .eq.
'GGA' .or. functl .eq.
'gga')
then
328 write(
message(1),
'(a,a)')
'atomxc: Unknown functional ', functl
334 call xc_f03_func_init(x_conf, xc_gga_x_pbe, nspin)
335 call xc_f03_func_init(c_conf, xc_gga_c_pbe, nspin)
337 call xc_f03_func_init(x_conf, xc_lda_x, nspin)
338 if (author .eq.
'CA' .or. author .eq.
'ca' .or. &
339 author .eq.
'PZ' .or. author .eq.
'pz')
then
340 call xc_f03_func_init(c_conf, xc_lda_c_pz, nspin)
341 else if (author .eq.
'PW92' .or. author .eq.
'pw92')
then
342 call xc_f03_func_init(c_conf, xc_lda_c_pw, nspin)
344 write(
message(1),
'(a,a)')
'LDAXC: Unknown author ', author
370 in1 = max( 1, ir-nn ) - ir
371 in2 = min( nr, ir+nn ) - ir
379 if (jn /= 0) dgdm(in) = dgdm(in) +
m_one / (0 - jn)
385 if (jn /= in .and. jn /= 0) f1 = f1 * (0 - jn)
386 if (jn /= in) f2 = f2 * (in - jn)
395 drdm = drdm + rmesh(ir+in) * dgdm(in)
399 dvol = 4 *
m_pi * rmesh(ir)**2 * drdm
402 dvol = safe_tol(dvol, dvmin)
403 if (ir == 1 .or. ir == nr) dvol = dvol / 2
404 if (gga) aux(ir) = dvol
410 dgidfj(in) = dgdm(in) / drdm
422 gd(3,is) = gd(3,is) + dgidfj(in) * dens(ir+in,is)
431 sigma(1) = gd(3, 1)*gd(3, 1)
433 sigma(2) = gd(3, 1) * gd(3, 2)
434 sigma(3) = gd(3, 2) * gd(3, 2)
440 call xc_f03_gga_exc_vxc(x_conf, int(1, c_size_t), d, sigma, epsx, dexdd, vxsigma)
441 call xc_f03_gga_exc_vxc(c_conf, int(1, c_size_t), d, sigma, epsc, decdd, vcsigma)
443 call xc_f03_lda_exc_vxc(x_conf, int(1, c_size_t), d, epsx, dexdd)
444 call xc_f03_lda_exc_vxc(c_conf, int(1, c_size_t), d, epsc, decdd)
451 dexdgd(:, 1) =
m_two * vxsigma(1) * gd(:, 1)
452 decdgd(:, 1) =
m_two * vcsigma(1) * gd(:, 1)
454 dexdgd(:, 1) = dexdgd(:, 1) + vxsigma(2)*gd(:, 2)
455 dexdgd(:, 2) =
m_two * vxsigma(3) * gd(:, 2) + vxsigma(2) * gd(:, 1)
456 decdgd(:, 1) = decdgd(:, 1) + vcsigma(2)*gd(:, 2)
457 decdgd(:, 2) =
m_two * vcsigma(3) * gd(:, 2) + vcsigma(2) * gd(:, 1)
464 ex = ex + dvol * d(is) * epsx(1)
465 ec = ec + dvol * d(is) * epsc(1)
466 dx = dx + dvol * d(is) * (epsx(1) - dexdd(is))
467 dc = dc + dvol * d(is) * (epsc(1) - decdd(is))
469 v_xc(ir,is) = v_xc(ir,is) + dvol * (dexdd(is) + decdd(is))
471 dx = dx - dvol * dens(ir+in,is) * dexdgd(3,is) * dgidfj(in)
472 dc = dc - dvol * dens(ir+in,is) * decdgd(3,is) * dgidfj(in)
473 v_xc(ir+in,is) = v_xc(ir+in,is) + dvol * &
474 (dexdgd(3,is) + decdgd(3,is)) * dgidfj(in)
477 v_xc(ir,is) = dexdd(is) + decdd(is)
491 v_xc(ir,is) = v_xc(ir,is) / dvol
506 v_xc(ir,is) = v_xc(ir,is) / eunit
510 call xc_f03_func_end(x_conf)
511 call xc_f03_func_end(c_conf)
536 subroutine vhrtre(rho, v, r, drdi, srdrdi, nr, a)
537 real(real64),
intent(in) :: rho(:)
538 real(real64),
intent(out) :: v(:)
539 real(real64),
intent(in) :: r(:), drdi(:), srdrdi(:)
540 integer,
intent(in) :: nr
541 real(real64),
intent(in) :: a
543 real(real64) :: x,y,q,a2by4,ybyq,qbyy,qpartc,v0,qt,dz,t,beta,dv
544 integer :: nrm1,nrm2,ir
551 ybyq =
m_one - a*a/48.0_real64
566 dz = drdi(ir)*rho(ir)
572 dz = drdi(ir)*rho(ir)
576 dz = drdi(nr)*rho(nr)
578 qt = (qt + qt + dz)*
m_half
579 v0 = (v0 + v0 + dz/r(nr))
590 beta = drdi(ir)*t*rho(ir)
608 q = (y - beta / 12.0_real64) * qbyy
609 v(ir) =
m_two * t * q
616 x = x + a2by4 * q - beta
618 t = srdrdi(ir) / r(ir)
619 beta = t * drdi(ir) * rho(ir)
620 q = (y - beta / 12.0_real64) * qbyy
621 v(ir) =
m_two * t * q
633 qpartc = r(nr)*v(nr)/
m_two
674 subroutine egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
676 real(real64),
intent(in) :: h(:), s(:)
677 integer,
intent(in) :: n
678 real(real64),
intent(inout) :: e
679 real(real64),
intent(out) :: g(:), y(:)
680 integer,
intent(in) :: l
681 real(real64),
intent(in) :: z, a, b, rmax
682 integer,
intent(in) :: nprin, nnode
683 real(real64),
intent(in) :: dr
684 integer,
intent(out) :: ierr
686 integer :: i,ncor,n1,n2,niter,nt
687 real(real64) :: e1, e2, del, de, et, t
688 real(real64),
parameter :: tol = 1.0e-5_real64
736 if (.not. (et <= e1 .or. et >= e2 .or. nt < nnode-1 .or. nt > nnode))
then
738 if (abs(de) < tol)
then
743 call yofe(e,de,dr,rmax,h,s,y,n,l,ncor,nt,z,a,b)
766 if (n2 >= nnode) cycle
786 if (n1 < nnode) cycle
805 g(i) = y(i) / (
m_one - t / 12._real64)
807 call nrmlzg(namespace, g, s, n)
814 subroutine yofe(e,de,dr,rmax,h,s,y,nmax,l,ncor,nnode,z,a,b)
835 real(real64),
intent(inout) :: e
836 real(real64),
intent(out) :: de
837 real(real64),
intent(in) :: dr, rmax
838 real(real64),
intent(in) :: h(:), s(:)
839 real(real64),
intent(out) :: y(:)
840 integer,
intent(in) :: nmax, l
841 integer,
intent(in) :: ncor
842 integer,
intent(out) :: nnode
843 real(real64),
intent(in) :: z, a, b
845 real(real64) :: y2, g, gsg, x, gin, gsgin, xin
846 integer :: n,knk,nndin,i
847 real(real64) :: zdr, yn, ratio, t
854 if (h(n)-e*s(n) <
m_one)
then
861 call bcorgn(e,h,s,l,zdr,y2)
870 call numout(e,h,s,y,ncor,knk,nnode,y2,g,gsg,x)
884 if (.not. (n < nmax .or. abs(dr) > 1.e3_real64))
then
885 call bcrmax(e,dr,rmax,h,s,n,yn,a,b)
899 call numin(e,h,s,y,n,nndin,yn,gin,gsgin,xin,knk)
913 gsg = gsg + gsgin * ratio * ratio
914 t = h(knk) - e * s(knk)
915 de = g * (x + xin + t * g) / gsg
916 nnode = nnode + nndin
917 if (de <
m_zero) nnode = nnode + 1
926 subroutine nrmlzg(namespace, g, s, n)
949 real(real64),
intent(inout) :: g(:)
950 real(real64),
intent(in) :: s(:)
951 integer,
intent(in) :: n
954 real(real64) :: norm, srnrm
963 if (mod(n,2) /= 1)
then
964 write(
message(1),
'(a,i6)')
' nrmlzg: n should be odd. n =', n
971 norm=norm + g(i)*s(i)*g(i)
976 norm=norm + g(i)*s(i)*g(i)
981 norm = norm + g(i) * s(i) * g(i)
993 subroutine bcorgn(e,h,s,l,zdr,y2)
994 real(real64),
intent(in) :: e
995 real(real64),
intent(in) :: h(:), s(:)
996 integer,
intent(in) :: l
997 real(real64),
intent(in) :: zdr
998 real(real64),
intent(out) :: y2
1000 real(real64) :: t2, t3, d2, c0, c1, c2
1009 t2 = h(2) - e * s(2)
1010 d2 = -((24._real64 + 10.0_real64 * t2) / (12._real64 - t2))
1026 c0=c0/(
m_one-0.75_real64*zdr)
1032 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1034 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1035 d2=(d2-c1)/(
m_one-c2)
1042 c0=c0/(
m_one-0.75_real64*zdr)
1047 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1049 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1050 d2=(d2-c1)/(
m_one-c2)
1060 subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
1061 real(real64),
intent(in) :: e, dr, rmax
1062 real(real64),
intent(in) :: h(:), s(:)
1063 integer,
intent(in) :: n
1064 real(real64),
intent(out) :: yn
1065 real(real64),
intent(in) :: a, b
1067 real(real64) :: tnm1, tn, tnp1, beta, dg, c1, c2, c3, dn
1071 tnm1=h(n-1)-e*s(n-1)
1073 tnp1=h(n+1)-e*s(n+1)
1077 c2=24._real64*dg/(12._real64-tn)
1078 dn=-((24._real64+10._real64*tn)/(12._real64-tn))
1080 c1= (
m_one-tnm1/6.0_real64)/(
m_one-tnm1/12._real64)
1081 c3=-((
m_one-tnp1/6.0_real64)/(
m_one-tnp1/12._real64))
1082 yn=-((
m_one-c1/c3)/(dn-c2/c3))
1088 subroutine numin(e,h,s,y,n,nnode,yn,g,gsg,x,knk)
1089 real(real64),
intent(in) :: e, h(:), s(:)
1090 real(real64),
intent(inout) :: y(:)
1091 integer,
intent(in) :: n
1092 integer,
intent(out) :: nnode
1093 real(real64),
intent(in) :: yn
1094 real(real64),
intent(out) :: g, gsg, x
1095 integer,
intent(in) :: knk
1104 g=y(n)/(
m_one-t/12._real64)
1109 g=y(i)/(
m_one-t/12._real64)
1123 do i = n - 2, knk, -1
1130 g = y(i) / (
m_one - t / 12._real64)
1131 gsg = gsg + g * s(i) * g
1145 end subroutine numin
1148 subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
1149 real(real64),
intent(in) :: e, h(:), s(:)
1150 real(real64),
intent(inout) :: y(:)
1151 integer,
intent(in) :: ncor
1152 integer,
intent(inout) :: knk
1153 integer,
intent(out) :: nnode
1154 real(real64),
intent(in) :: y2
1155 real(real64),
intent(out) :: g, gsg, x
1158 real(real64) :: t, xl
1165 g=y(2)/(
m_one-t/12._real64)
1169 g=y(3)/(
m_one-t/12._real64)
1190 g = y(i) / (
m_one - t / 12._real64)
1191 gsg = gsg + g * s(i) * g
1192 if (.not. (nnode < ncor .or. xl*x >
m_zero))
then
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, public write_valconf(c, s)
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
subroutine, public valconf_copy(cout, cin)
subroutine, public read_valconf(namespace, s, c)
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)