25 use,
intrinsic :: iso_fortran_env
46 valconf_unpolarized_to_polarized, &
49 character(len=1),
parameter :: &
50 spec_notation(0:3) = (/
's',
'p',
'd',
'f' /)
52 integer,
parameter :: VALCONF_STRING_LENGTH = 80
57 character(len=3) :: symbol =
""
62 real(real64) :: occ(15,2) =
m_zero
63 real(real64) :: j(15) =
m_zero
70 subroutine valconf_copy(cout, cin)
71 type(valconf_t),
intent(out) :: cout
72 type(valconf_t),
intent(in) :: cin
74 push_sub(valconf_copy)
77 cout%symbol = cin%symbol
86 end subroutine valconf_copy
90 subroutine write_valconf(c, s)
91 type(valconf_t),
intent(in) :: c
92 character(len=VALCONF_STRING_LENGTH),
intent(out) :: s
96 push_sub(write_valconf)
98 write(s,
'(i2,1x,a2,i1,1x,i1,a1,6(i1,a1,f6.3,a1))') c%z, c%symbol, c%type, c%p,
':',&
99 (c%n(j),spec_notation(c%l(j)),c%occ(j,1),
',',j=1,c%p)
101 pop_sub(write_valconf)
102 end subroutine write_valconf
106 subroutine read_valconf(namespace, s, c)
107 type(namespace_t),
intent(in) :: namespace
108 character(len=VALCONF_STRING_LENGTH),
intent(in) :: s
109 type(valconf_t),
intent(out) :: c
112 character(len=1) :: lvalues(1:6)
114 push_sub(read_valconf)
116 read (s,
'(i2,1x,a2,i1,1x,i1,1x,6(i1,a1,f6.3,1x))') c%z, c%symbol, c%type, c%p,&
117 (c%n(j),lvalues(j),c%occ(j,1),j=1,c%p)
119 select case (lvalues(j))
136 pop_sub(read_valconf)
137 end subroutine read_valconf
140 subroutine valconf_unpolarized_to_polarized(conf)
141 type(valconf_t),
intent(inout) :: conf
146 push_sub(valconf_unpolarized_to_polarized)
151 conf%occ(ii, 1) = min(xx, real(2*ll+1, real64) )
152 conf%occ(ii, 2) = xx - conf%occ(ii,1)
155 pop_sub(valconf_unpolarized_to_polarized)
156 end subroutine valconf_unpolarized_to_polarized
160 subroutine atomhxc(namespace, functl, g, nspin, dens, v, extra)
161 type(namespace_t),
intent(in) :: namespace
162 character(len=*),
intent(in) :: functl
163 type(logrid_t),
intent(in) :: g
164 integer,
intent(in) :: nspin
165 real(real64),
intent(in) :: dens(g%nrval, nspin)
166 real(real64),
intent(out) :: v(g%nrval, nspin)
167 real(real64),
intent(in),
optional :: extra(g%nrval)
169 character(len=5) :: xcfunc, xcauth
171 real(real64),
allocatable :: xc(:, :), ve(:, :), rho(:, :)
172 real(real64) :: r2, ex, ec, dx, dc
176 safe_allocate( ve(1:g%nrval, 1:nspin))
177 safe_allocate( xc(1:g%nrval, 1:nspin))
178 safe_allocate(rho(1:g%nrval, 1:nspin))
185 rho(:, 1) = rho(:, 1) + dens(:, is)
189 call vhrtre(rho(:, 1), ve(:, is), g%rofi, g%drdi, g%s, g%nrval, g%a)
191 ve(1, 1:nspin) = ve(2, 1:nspin)
201 message(1) =
'Internal Error in atomhxc: unknown functl'
207 if (
present(extra)) rho(:, is) = rho(:, is) + extra(:)/nspin
208 rho(2:g%nrval, is) = rho(2:g%nrval, is)/(
m_four*
m_pi*g%r2ofi(2:g%nrval))
210 r2 = g%rofi(2)/(g%rofi(3)-g%rofi(2))
211 rho(1, 1:nspin) = rho(2, 1:nspin) - (rho(3, 1:nspin)-rho(2, 1:nspin))*r2
219 call atomxc(namespace, xcfunc, xcauth, g%nrval, g%nrval, g%rofi, nspin, rho, ex, ec, dx, dc, xc)
222 safe_deallocate_a(ve)
223 safe_deallocate_a(xc)
224 safe_deallocate_a(rho)
227 end subroutine atomhxc
269 subroutine atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
270 type(namespace_t),
intent(in) :: namespace
271 character(len=*),
intent(in) :: FUNCTL, AUTHOR
272 integer,
intent(in) :: NR, MAXR
273 real(real64),
intent(in) :: RMESH(MAXR)
274 integer,
intent(in) :: NSPIN
275 real(real64),
intent(in) :: DENS(MAXR,NSPIN)
276 real(real64),
intent(out) :: DC, DX, EC, EX, V_XC(MAXR,NSPIN)
284 integer,
parameter :: NN = 5
292 real(real64),
parameter :: EUNIT =
m_half
298 real(real64),
parameter :: DVMIN = 1.0e-12_real64
305 integer :: IN, IN1, IN2, IR, IS, JN
306 real(real64) :: AUX(MAXR), D(NSPIN), &
307 DGDM(-NN:NN), DGIDFJ(-NN:NN), DRDM, DVOL, &
308 F1, F2, GD(3,NSPIN), sigma(3), vxsigma(3), vcsigma(3), &
310 DEXDD(NSPIN), DEXDGD(3,NSPIN), DECDD(NSPIN), DECDGD(3,NSPIN)
311 type(xc_f03_func_t) :: x_conf, c_conf
316 assert(nspin == 1 .or. nspin == 2)
319 if (functl .eq.
'LDA' .or. functl .eq.
'lda' .or. &
320 functl .eq.
'LSD' .or. functl .eq.
'lsd')
then
322 else if (functl .eq.
'GGA' .or. functl .eq.
'gga')
then
326 write(
message(1),
'(a,a)')
'atomxc: Unknown functional ', functl
332 call xc_f03_func_init(x_conf, xc_gga_x_pbe, nspin)
333 call xc_f03_func_init(c_conf, xc_gga_c_pbe, nspin)
335 call xc_f03_func_init(x_conf, xc_lda_x, nspin)
336 if (author .eq.
'CA' .or. author .eq.
'ca' .or. &
337 author .eq.
'PZ' .or. author .eq.
'pz')
then
338 call xc_f03_func_init(c_conf, xc_lda_c_pz, nspin)
339 else if (author .eq.
'PW92' .or. author .eq.
'pw92')
then
340 call xc_f03_func_init(c_conf, xc_lda_c_pw, nspin)
342 write(
message(1),
'(a,a)')
'LDAXC: Unknown author ', author
368 in1 = max( 1, ir-nn ) - ir
369 in2 = min( nr, ir+nn ) - ir
377 if (jn /= 0) dgdm(in) = dgdm(in) +
m_one / (0 - jn)
383 if (jn /= in .and. jn /= 0) f1 = f1 * (0 - jn)
384 if (jn /= in) f2 = f2 * (in - jn)
393 drdm = drdm + rmesh(ir+in) * dgdm(in)
397 dvol = 4 *
m_pi * rmesh(ir)**2 * drdm
400 dvol = safe_tol(dvol, dvmin)
401 if (ir == 1 .or. ir == nr) dvol = dvol / 2
402 if (gga) aux(ir) = dvol
408 dgidfj(in) = dgdm(in) / drdm
420 gd(3,is) = gd(3,is) + dgidfj(in) * dens(ir+in,is)
429 sigma(1) = gd(3, 1)*gd(3, 1)
431 sigma(2) = gd(3, 1) * gd(3, 2)
432 sigma(3) = gd(3, 2) * gd(3, 2)
438 call xc_f03_gga_exc_vxc(x_conf, int(1, c_size_t), d, sigma, epsx, dexdd, vxsigma)
439 call xc_f03_gga_exc_vxc(c_conf, int(1, c_size_t), d, sigma, epsc, decdd, vcsigma)
441 call xc_f03_lda_exc_vxc(x_conf, int(1, c_size_t), d, epsx, dexdd)
442 call xc_f03_lda_exc_vxc(c_conf, int(1, c_size_t), d, epsc, decdd)
449 dexdgd(:, 1) =
m_two * vxsigma(1) * gd(:, 1)
450 decdgd(:, 1) =
m_two * vcsigma(1) * gd(:, 1)
452 dexdgd(:, 1) = dexdgd(:, 1) + vxsigma(2)*gd(:, 2)
453 dexdgd(:, 2) =
m_two * vxsigma(3) * gd(:, 2) + vxsigma(2) * gd(:, 1)
454 decdgd(:, 1) = decdgd(:, 1) + vcsigma(2)*gd(:, 2)
455 decdgd(:, 2) =
m_two * vcsigma(3) * gd(:, 2) + vcsigma(2) * gd(:, 1)
462 ex = ex + dvol * d(is) * epsx(1)
463 ec = ec + dvol * d(is) * epsc(1)
464 dx = dx + dvol * d(is) * (epsx(1) - dexdd(is))
465 dc = dc + dvol * d(is) * (epsc(1) - decdd(is))
467 v_xc(ir,is) = v_xc(ir,is) + dvol * (dexdd(is) + decdd(is))
469 dx = dx - dvol * dens(ir+in,is) * dexdgd(3,is) * dgidfj(in)
470 dc = dc - dvol * dens(ir+in,is) * decdgd(3,is) * dgidfj(in)
471 v_xc(ir+in,is) = v_xc(ir+in,is) + dvol * &
472 (dexdgd(3,is) + decdgd(3,is)) * dgidfj(in)
475 v_xc(ir,is) = dexdd(is) + decdd(is)
489 v_xc(ir,is) = v_xc(ir,is) / dvol
504 v_xc(ir,is) = v_xc(ir,is) / eunit
508 call xc_f03_func_end(x_conf)
509 call xc_f03_func_end(c_conf)
512 end subroutine atomxc
534 subroutine vhrtre(rho, v, r, drdi, srdrdi, nr, a)
535 real(real64),
intent(in) :: rho(:)
536 real(real64),
intent(out) :: v(:)
537 real(real64),
intent(in) :: r(:), drdi(:), srdrdi(:)
538 integer,
intent(in) :: nr
539 real(real64),
intent(in) :: a
541 real(real64) :: x,y,q,a2by4,ybyq,qbyy,qpartc,v0,qt,dz,t,beta,dv
542 integer :: nrm1,nrm2,ir
549 ybyq =
m_one - a*a/48.0_real64
564 dz = drdi(ir)*rho(ir)
570 dz = drdi(ir)*rho(ir)
574 dz = drdi(nr)*rho(nr)
576 qt = (qt + qt + dz)*
m_half
577 v0 = (v0 + v0 + dz/r(nr))
588 beta = drdi(ir)*t*rho(ir)
606 q = (y - beta / 12.0_real64) * qbyy
607 v(ir) =
m_two * t * q
614 x = x + a2by4 * q - beta
616 t = srdrdi(ir) / r(ir)
617 beta = t * drdi(ir) * rho(ir)
618 q = (y - beta / 12.0_real64) * qbyy
619 v(ir) =
m_two * t * q
631 qpartc = r(nr)*v(nr)/
m_two
646 end subroutine vhrtre
672 subroutine egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
673 type(namespace_t),
intent(in) :: namespace
674 real(real64),
intent(in) :: h(:), s(:)
675 integer,
intent(in) :: n
676 real(real64),
intent(inout) :: e
677 real(real64),
intent(out) :: g(:), y(:)
678 integer,
intent(in) :: l
679 real(real64),
intent(in) :: z, a, b, rmax
680 integer,
intent(in) :: nprin, nnode
681 real(real64),
intent(in) :: dr
682 integer,
intent(out) :: ierr
684 integer :: i,ncor,n1,n2,niter,nt
685 real(real64) :: e1, e2, del, de, et, t
686 real(real64),
parameter :: tol = 1.0e-5_real64
734 if (.not. (et <= e1 .or. et >= e2 .or. nt < nnode-1 .or. nt > nnode))
then
736 if (abs(de) < tol)
then
741 call yofe(e,de,dr,rmax,h,s,y,n,l,ncor,nt,z,a,b)
764 if (n2 >= nnode) cycle
784 if (n1 < nnode) cycle
803 g(i) = y(i) / (
m_one - t / 12._real64)
805 call nrmlzg(namespace, g, s, n)
812 subroutine yofe(e,de,dr,rmax,h,s,y,nmax,l,ncor,nnode,z,a,b)
833 real(real64),
intent(inout) :: e
834 real(real64),
intent(out) :: de
835 real(real64),
intent(in) :: dr, rmax
836 real(real64),
intent(in) :: h(:), s(:)
837 real(real64),
intent(out) :: y(:)
838 integer,
intent(in) :: nmax, l
839 integer,
intent(in) :: ncor
840 integer,
intent(out) :: nnode
841 real(real64),
intent(in) :: z, a, b
843 real(real64) :: y2, g, gsg, x, gin, gsgin, xin
844 integer :: n,knk,nndin,i
845 real(real64) :: zdr, yn, ratio, t
852 if (h(n)-e*s(n) <
m_one)
then
859 call bcorgn(e,h,s,l,zdr,y2)
868 call numout(e,h,s,y,ncor,knk,nnode,y2,g,gsg,x)
882 if (.not. (n < nmax .or. abs(dr) > 1.e3_real64))
then
883 call bcrmax(e,dr,rmax,h,s,n,yn,a,b)
897 call numin(e,h,s,y,n,nndin,yn,gin,gsgin,xin,knk)
911 gsg = gsg + gsgin * ratio * ratio
912 t = h(knk) - e * s(knk)
913 de = g * (x + xin + t * g) / gsg
914 nnode = nnode + nndin
915 if (de <
m_zero) nnode = nnode + 1
924 subroutine nrmlzg(namespace, g, s, n)
946 type(namespace_t),
intent(in) :: namespace
947 real(real64),
intent(inout) :: g(:)
948 real(real64),
intent(in) :: s(:)
949 integer,
intent(in) :: n
952 real(real64) :: norm, srnrm
961 if (mod(n,2) /= 1)
then
962 write(
message(1),
'(a,i6)')
' nrmlzg: n should be odd. n =', n
969 norm=norm + g(i)*s(i)*g(i)
974 norm=norm + g(i)*s(i)*g(i)
979 norm = norm + g(i) * s(i) * g(i)
988 end subroutine nrmlzg
991 subroutine bcorgn(e,h,s,l,zdr,y2)
992 real(real64),
intent(in) :: e
993 real(real64),
intent(in) :: h(:), s(:)
994 integer,
intent(in) :: l
995 real(real64),
intent(in) :: zdr
996 real(real64),
intent(out) :: y2
998 real(real64) :: t2, t3, d2, c0, c1, c2
1007 t2 = h(2) - e * s(2)
1008 d2 = -((24._real64 + 10.0_real64 * t2) / (12._real64 - t2))
1024 c0=c0/(
m_one-0.75_real64*zdr)
1030 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1032 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1033 d2=(d2-c1)/(
m_one-c2)
1040 c0=c0/(
m_one-0.75_real64*zdr)
1045 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1047 c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
1048 d2=(d2-c1)/(
m_one-c2)
1052 end subroutine bcorgn
1058 subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
1059 real(real64),
intent(in) :: e, dr, rmax
1060 real(real64),
intent(in) :: h(:), s(:)
1061 integer,
intent(in) :: n
1062 real(real64),
intent(out) :: yn
1063 real(real64),
intent(in) :: a, b
1065 real(real64) :: tnm1, tn, tnp1, beta, dg, c1, c2, c3, dn
1069 tnm1=h(n-1)-e*s(n-1)
1071 tnp1=h(n+1)-e*s(n+1)
1075 c2=24._real64*dg/(12._real64-tn)
1076 dn=-((24._real64+10._real64*tn)/(12._real64-tn))
1078 c1= (
m_one-tnm1/6.0_real64)/(
m_one-tnm1/12._real64)
1079 c3=-((
m_one-tnp1/6.0_real64)/(
m_one-tnp1/12._real64))
1080 yn=-((
m_one-c1/c3)/(dn-c2/c3))
1083 end subroutine bcrmax
1086 subroutine numin(e,h,s,y,n,nnode,yn,g,gsg,x,knk)
1087 real(real64),
intent(in) :: e, h(:), s(:)
1088 real(real64),
intent(inout) :: y(:)
1089 integer,
intent(in) :: n
1090 integer,
intent(out) :: nnode
1091 real(real64),
intent(in) :: yn
1092 real(real64),
intent(out) :: g, gsg, x
1093 integer,
intent(in) :: knk
1102 g=y(n)/(
m_one-t/12._real64)
1107 g=y(i)/(
m_one-t/12._real64)
1121 do i = n - 2, knk, -1
1128 g = y(i) / (
m_one - t / 12._real64)
1129 gsg = gsg + g * s(i) * g
1143 end subroutine numin
1146 subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
1147 real(real64),
intent(in) :: e, h(:), s(:)
1148 real(real64),
intent(inout) :: y(:)
1149 integer,
intent(in) :: ncor
1150 integer,
intent(inout) :: knk
1151 integer,
intent(out) :: nnode
1152 real(real64),
intent(in) :: y2
1153 real(real64),
intent(out) :: g, gsg, x
1156 real(real64) :: t, xl
1163 g=y(2)/(
m_one-t/12._real64)
1167 g=y(3)/(
m_one-t/12._real64)
1188 g = y(i) / (
m_one - t / 12._real64)
1189 gsg = gsg + g * s(i) * g
1190 if (.not. (nnode < ncor .or. xl*x >
m_zero))
then
1201 end subroutine numout
1203end module atomic_oct_m
double sqrt(double __x) __attribute__((__nothrow__
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)