24  use, 
intrinsic :: iso_fortran_env
 
   48  character(len=1), 
parameter :: &
 
   49    spec_notation(0:3) = (/ 
's', 
'p', 
'd', 
'f' /)
 
   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
 
   90    type(valconf_t), 
intent(in) :: c
 
   91    character(len=VALCONF_STRING_LENGTH), 
intent(out) :: s
 
   97    write(s,
'(i2,1x,a2,i1,1x,i1,a1,6(i1,a1,f6.3,a1))') c%z, c%symbol, c%type, c%p, 
':',&
 
   98      (c%n(j),spec_notation(c%l(j)),c%occ(j,1),
',',j=1,c%p)
 
  106    type(namespace_t),        
intent(in)    :: namespace
 
  107    character(len=VALCONF_STRING_LENGTH), 
intent(in) :: s
 
  108    type(valconf_t), 
intent(out) :: c
 
  111    character(len=1) :: lvalues(1:6)
 
  115    read (s,
'(i2,1x,a2,i1,1x,i1,1x,6(i1,a1,f6.3,1x))') c%z, c%symbol, c%type, c%p,&
 
  116      (c%n(j),lvalues(j),c%occ(j,1),j=1,c%p)
 
  118      select case (lvalues(j))
 
  140    type(valconf_t), 
intent(inout) :: conf
 
  150      conf%occ(ii, 1) = min(xx, real(2*ll+1, real64) )
 
  151      conf%occ(ii, 2) = xx - conf%occ(ii,1)
 
  159  subroutine atomhxc(namespace, functl, g, nspin, dens, v, extra)
 
  161    character(len=*),  
intent(in)  :: functl
 
  163    integer,           
intent(in)  :: nspin
 
  164    real(real64),      
intent(in)  :: dens(g%nrval, nspin)
 
  165    real(real64),      
intent(out) :: v(g%nrval, nspin)
 
  166    real(real64),      
intent(in), 
optional :: extra(g%nrval)
 
  168    character(len=5) :: xcfunc, xcauth
 
  170    real(real64), 
allocatable :: xc(:, :), ve(:, :), rho(:, :)
 
  171    real(real64) :: r2, ex, ec, dx, dc
 
  175    safe_allocate( ve(1:g%nrval, 1:nspin))
 
  176    safe_allocate( xc(1:g%nrval, 1:nspin))
 
  177    safe_allocate(rho(1:g%nrval, 1:nspin))
 
  184      rho(:, 1) = rho(:, 1) + dens(:, is)
 
  188      call vhrtre(rho(:, 1), ve(:, is), g%rofi, g%drdi, g%s, g%nrval, g%a)
 
  190    ve(1, 1:nspin) = ve(2, 1:nspin)
 
  200      message(1) = 
'Internal Error in atomhxc: unknown functl' 
  206      if (
present(extra)) rho(:, is) = rho(:, is) + extra(:)/nspin
 
  207      rho(2:g%nrval, is) = rho(2:g%nrval, is)/(
m_four*
m_pi*g%r2ofi(2:g%nrval))
 
  209    r2 = g%rofi(2)/(g%rofi(3)-g%rofi(2))
 
  210    rho(1, 1:nspin) = rho(2, 1:nspin) - (rho(3, 1:nspin)-rho(2, 1:nspin))*r2
 
  218    call atomxc(namespace, xcfunc, xcauth, g%nrval, g%nrval, g%rofi, nspin, rho, ex, ec, dx, dc, xc)
 
  221    safe_deallocate_a(ve)
 
  222    safe_deallocate_a(xc)
 
  223    safe_deallocate_a(rho)
 
  268  subroutine atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
 
  270    character(len=*), 
intent(in) :: functl, author
 
  271    integer, 
intent(in) :: nr, maxr
 
  272    real(real64), 
intent(in) :: rmesh(maxr)
 
  273    integer, 
intent(in) :: nspin
 
  274    real(real64), 
intent(in) :: dens(maxr,nspin)
 
  275    real(real64), 
intent(out) :: dc, dx, ec, ex, v_xc(maxr,nspin)
 
  283    integer, 
parameter :: nn = 5
 
  291    real(real64), 
parameter :: eunit = 
m_half 
  297    real(real64), 
parameter :: dvmin = 1.0e-12_real64
 
  304    integer :: in, in1, in2, ir, is, jn
 
  305    real(real64) :: aux(maxr), d(nspin),                              &
 
  306      dgdm(-nn:nn), dgidfj(-nn:nn), drdm, dvol,                      &
 
  307      f1, f2, gd(3,nspin), sigma(3), vxsigma(3), vcsigma(3),   &
 
  309      dexdd(nspin), dexdgd(3,nspin), decdd(nspin), decdgd(3,nspin)
 
  310    type(xc_f03_func_t) :: x_conf, c_conf
 
  315    assert(nspin == 1 .or. nspin == 2)
 
  318    if (functl .eq. 
'LDA' .or. functl .eq. 
'lda' .or.                              &
 
  319      functl .eq. 
'LSD' .or. functl .eq. 
'lsd') 
then 
  321    else if (functl .eq. 
'GGA' .or. functl .eq. 
'gga') 
then 
  325      write(
message(1),
'(a,a)') 
'atomxc: Unknown functional ', functl
 
  331      call xc_f03_func_init(x_conf, xc_gga_x_pbe, nspin)
 
  332      call xc_f03_func_init(c_conf, xc_gga_c_pbe, nspin)
 
  334      call xc_f03_func_init(x_conf, xc_lda_x, nspin)
 
  335      if (author .eq. 
'CA' .or. author .eq. 
'ca' .or.                                &
 
  336        author .eq. 
'PZ' .or. author .eq. 
'pz') 
then 
  337        call xc_f03_func_init(c_conf, xc_lda_c_pz, nspin)
 
  338      else if (author .eq. 
'PW92' .or. author .eq. 
'pw92') 
then 
  339        call xc_f03_func_init(c_conf, xc_lda_c_pw, nspin)
 
  341        write(
message(1),
'(a,a)') 
'LDAXC: Unknown author ', author
 
  367      in1 = max(  1, ir-nn ) - ir
 
  368      in2 = min( nr, ir+nn ) - ir
 
  376            if (jn /= 0) dgdm(in) = dgdm(in) + 
m_one / (0 - jn)
 
  382            if (jn /= in .and. jn /= 0) f1 = f1 * (0  - jn)
 
  383            if (jn /= in)               f2 = f2 * (in - jn)
 
  392        drdm = drdm + rmesh(ir+in) * dgdm(in)
 
  396      dvol = 4 * 
m_pi * rmesh(ir)**2 * drdm
 
  399      dvol = safe_tol(dvol, dvmin)
 
  400      if (ir == 1 .or. ir == nr) dvol = dvol / 2
 
  401      if (gga) aux(ir) = dvol
 
  407          dgidfj(in) = dgdm(in) / drdm
 
  419            gd(3,is) = gd(3,is) + dgidfj(in) * dens(ir+in,is)
 
  428      sigma(1) = gd(3, 1)*gd(3, 1)
 
  430        sigma(2) = gd(3, 1) * gd(3, 2)
 
  431        sigma(3) = gd(3, 2) * gd(3, 2)
 
  437        call xc_f03_gga_exc_vxc(x_conf, int(1, c_size_t), d, sigma, epsx, dexdd, vxsigma)
 
  438        call xc_f03_gga_exc_vxc(c_conf, int(1, c_size_t), d, sigma, epsc, decdd, vcsigma)
 
  440        call xc_f03_lda_exc_vxc(x_conf, int(1, c_size_t), d, epsx, dexdd)
 
  441        call xc_f03_lda_exc_vxc(c_conf, int(1, c_size_t), d, epsc, decdd)
 
  448        dexdgd(:, 1) = 
m_two * vxsigma(1) * gd(:, 1)
 
  449        decdgd(:, 1) = 
m_two * vcsigma(1) * gd(:, 1)
 
  451          dexdgd(:, 1) = dexdgd(:, 1) + vxsigma(2)*gd(:, 2)
 
  452          dexdgd(:, 2) = 
m_two * vxsigma(3) * gd(:, 2) + vxsigma(2) * gd(:, 1)
 
  453          decdgd(:, 1) = decdgd(:, 1) + vcsigma(2)*gd(:, 2)
 
  454          decdgd(:, 2) = 
m_two * vcsigma(3) * gd(:, 2) + vcsigma(2) * gd(:, 1)
 
  461        ex = ex + dvol * d(is) * epsx(1)
 
  462        ec = ec + dvol * d(is) * epsc(1)
 
  463        dx = dx + dvol * d(is) * (epsx(1) - dexdd(is))
 
  464        dc = dc + dvol * d(is) * (epsc(1) - decdd(is))
 
  466          v_xc(ir,is) = v_xc(ir,is) + dvol * (dexdd(is) + decdd(is))
 
  468            dx = dx - dvol * dens(ir+in,is) * dexdgd(3,is) * dgidfj(in)
 
  469            dc = dc - dvol * dens(ir+in,is) * decdgd(3,is) * dgidfj(in)
 
  470            v_xc(ir+in,is) = v_xc(ir+in,is) + dvol *                         &
 
  471              (dexdgd(3,is) + decdgd(3,is)) * dgidfj(in)
 
  474          v_xc(ir,is) = dexdd(is) + decdd(is)
 
  488          v_xc(ir,is) = v_xc(ir,is) / dvol
 
  503        v_xc(ir,is) = v_xc(ir,is) / eunit
 
  507    call xc_f03_func_end(x_conf)
 
  508    call xc_f03_func_end(c_conf)
 
  533  subroutine vhrtre(rho, v, r, drdi, srdrdi, nr, a)
 
  534    real(real64),   
intent(in)  :: rho(:)
 
  535    real(real64),   
intent(out) :: v(:)
 
  536    real(real64),   
intent(in)  :: r(:), drdi(:), srdrdi(:)
 
  537    integer, 
intent(in)  :: nr
 
  538    real(real64),   
intent(in)  :: a
 
  540    real(real64) :: x,y,q,a2by4,ybyq,qbyy,qpartc,v0,qt,dz,t,beta,dv
 
  541    integer :: nrm1,nrm2,ir
 
  548    ybyq = 
m_one - a*a/48.0_real64
 
  563      dz = drdi(ir)*rho(ir)
 
  569      dz = drdi(ir)*rho(ir)
 
  573    dz = drdi(nr)*rho(nr)
 
  575    qt = (qt + qt + dz)*
m_half 
  576    v0 = (v0 + v0 + dz/r(nr))
 
  587    beta = drdi(ir)*t*rho(ir)
 
  605    q = (y - beta / 12.0_real64) * qbyy
 
  606    v(ir) = 
m_two * t * q
 
  613      x  = x + a2by4 * q - beta
 
  615      t  = srdrdi(ir) / r(ir)
 
  616      beta = t * drdi(ir) * rho(ir)
 
  617      q = (y - beta / 12.0_real64) * qbyy
 
  618      v(ir) = 
m_two * t * q
 
  630    qpartc = r(nr)*v(nr)/
m_two 
  671  subroutine egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
 
  673    real(real64),          
intent(in)    :: h(:), s(:)
 
  674    integer,           
intent(in)    :: n
 
  675    real(real64),          
intent(inout) :: e
 
  676    real(real64),          
intent(out)   :: g(:), y(:)
 
  677    integer,           
intent(in)    :: l
 
  678    real(real64),          
intent(in)    :: z, a, b, rmax
 
  679    integer,           
intent(in)    :: nprin, nnode
 
  680    real(real64),          
intent(in)    :: dr
 
  681    integer,           
intent(out)   :: ierr
 
  683    integer :: i,ncor,n1,n2,niter,nt
 
  684    real(real64)               :: e1, e2, del, de, et, t
 
  685    real(real64), 
parameter    :: tol = 1.0e-5_real64
 
  733      if (.not. (et <= e1 .or. et >= e2 .or. nt < nnode-1 .or. nt > nnode)) 
then 
  735        if (abs(de) < tol) 
then 
  740      call yofe(e,de,dr,rmax,h,s,y,n,l,ncor,nt,z,a,b)
 
  763        if (n2 >= nnode) cycle
 
  783      if (n1  <  nnode) cycle
 
  802        g(i) = y(i) / (
m_one - t / 12._real64)
 
  804      call nrmlzg(namespace, g, s, n)
 
  811  subroutine yofe(e,de,dr,rmax,h,s,y,nmax,l,ncor,nnode,z,a,b)
 
  832    real(real64), 
intent(inout) :: e
 
  833    real(real64), 
intent(out)   :: de
 
  834    real(real64), 
intent(in)    :: dr, rmax
 
  835    real(real64), 
intent(in)    :: h(:), s(:)
 
  836    real(real64), 
intent(out)   :: y(:)
 
  837    integer,  
intent(in)    :: nmax, l
 
  838    integer,  
intent(in)    :: ncor
 
  839    integer,  
intent(out)   :: nnode
 
  840    real(real64), 
intent(in)    :: z, a, b
 
  842    real(real64) :: y2, g, gsg, x, gin, gsgin, xin
 
  843    integer :: n,knk,nndin,i
 
  844    real(real64) :: zdr, yn, ratio, t
 
  851      if (h(n)-e*s(n) < 
m_one) 
then 
  858    call bcorgn(e,h,s,l,zdr,y2)
 
  867    call numout(e,h,s,y,ncor,knk,nnode,y2,g,gsg,x)
 
  881    if (.not. (n < nmax .or. abs(dr) > 1.e3_real64)) 
then 
  882      call bcrmax(e,dr,rmax,h,s,n,yn,a,b)
 
  896    call numin(e,h,s,y,n,nndin,yn,gin,gsgin,xin,knk)
 
  910    gsg = gsg + gsgin * ratio * ratio
 
  911    t = h(knk) - e * s(knk)
 
  912    de = g * (x + xin + t * g) / gsg
 
  913    nnode = nnode + nndin
 
  914    if (de < 
m_zero) nnode = nnode + 1
 
  923  subroutine nrmlzg(namespace, g, s, n)
 
  946    real(real64),          
intent(inout) :: g(:)
 
  947    real(real64),          
intent(in)    :: s(:)
 
  948    integer,           
intent(in)    :: n
 
  951    real(real64) :: norm, srnrm
 
  960    if (mod(n,2) /= 1) 
then 
  961      write(
message(1),
'(a,i6)') 
' nrmlzg: n should be odd. n =', n
 
  968      norm=norm + g(i)*s(i)*g(i)
 
  973      norm=norm + g(i)*s(i)*g(i)
 
  978      norm = norm + g(i) * s(i) * g(i)
 
  990  subroutine bcorgn(e,h,s,l,zdr,y2)
 
  991    real(real64), 
intent(in) :: e
 
  992    real(real64), 
intent(in) :: h(:), s(:)
 
  993    integer,     
intent(in) :: l
 
  994    real(real64), 
intent(in) :: zdr
 
  995    real(real64), 
intent(out) :: y2
 
  997    real(real64) :: t2, t3, d2, c0, c1, c2
 
 1006    t2 = h(2) - e * s(2)
 
 1007    d2 = -((24._real64 + 10.0_real64 * t2) / (12._real64 - t2))
 
 1023        c0=c0/(
m_one-0.75_real64*zdr)
 
 1029      c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
 
 1031      c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
 
 1032      d2=(d2-c1)/(
m_one-c2)
 
 1039        c0=c0/(
m_one-0.75_real64*zdr)
 
 1044      c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
 
 1046      c2=(-
m_half)*c0*(24._real64-t3)/(12._real64-t3)
 
 1047      d2=(d2-c1)/(
m_one-c2)
 
 1057  subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
 
 1058    real(real64), 
intent(in)  :: e, dr, rmax
 
 1059    real(real64), 
intent(in)  :: h(:), s(:)
 
 1060    integer,  
intent(in)  :: n
 
 1061    real(real64), 
intent(out) :: yn
 
 1062    real(real64), 
intent(in)  :: a, b
 
 1064    real(real64) :: tnm1, tn, tnp1, beta, dg, c1, c2, c3, dn
 
 1068    tnm1=h(n-1)-e*s(n-1)
 
 1070    tnp1=h(n+1)-e*s(n+1)
 
 1074    c2=24._real64*dg/(12._real64-tn)
 
 1075    dn=-((24._real64+10._real64*tn)/(12._real64-tn))
 
 1077    c1= (
m_one-tnm1/6.0_real64)/(
m_one-tnm1/12._real64)
 
 1078    c3=-((
m_one-tnp1/6.0_real64)/(
m_one-tnp1/12._real64))
 
 1079    yn=-((
m_one-c1/c3)/(dn-c2/c3))
 
 1085  subroutine numin(e,h,s,y,n,nnode,yn,g,gsg,x,knk)
 
 1086    real(real64), 
intent(in)    :: e, h(:), s(:)
 
 1087    real(real64), 
intent(inout) :: y(:)
 
 1088    integer,  
intent(in)    :: n
 
 1089    integer,  
intent(out)   :: nnode
 
 1090    real(real64), 
intent(in)    :: yn
 
 1091    real(real64), 
intent(out)   :: g, gsg, x
 
 1092    integer,  
intent(in)    :: knk
 
 1101    g=y(n)/(
m_one-t/12._real64)
 
 1106    g=y(i)/(
m_one-t/12._real64)
 
 1120    do i = n - 2, knk, -1
 
 1127      g = y(i) / (
m_one - t / 12._real64)
 
 1128      gsg = gsg + g * s(i) * g
 
 1142  end subroutine numin 
 1145  subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
 
 1146    real(real64), 
intent(in)    :: e, h(:), s(:)
 
 1147    real(real64), 
intent(inout) :: y(:)
 
 1148    integer,  
intent(in)    :: ncor
 
 1149    integer,  
intent(inout) :: knk
 
 1150    integer,  
intent(out)   :: nnode
 
 1151    real(real64), 
intent(in)    :: y2
 
 1152    real(real64), 
intent(out)   :: g, gsg, x
 
 1155    real(real64) :: t, xl
 
 1162    g=y(2)/(
m_one-t/12._real64)
 
 1166    g=y(3)/(
m_one-t/12._real64)
 
 1187      g = y(i) / (
m_one - t / 12._real64)
 
 1188      gsg = gsg + g * s(i) * g
 
 1189      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)