The Elk Code
atom.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: atom
8 ! !INTERFACE:
9 subroutine atom(sol,ptnucl,zn,nst,n,l,k,occ,xctype,xcgrad,nr,r,eval,rho,vr,rwf)
10 ! !USES:
11 use modxcifc
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! sol : speed of light in atomic units (in,real)
14 ! ptnucl : .true. if the nucleus is a point particle (in,logical)
15 ! zn : nuclear charge (in,real)
16 ! nst : number of states to solve for (in,integer)
17 ! n : priciple quantum number of each state (in,integer(nst))
18 ! l : quantum number l of each state (in,integer(nst))
19 ! k : quantum number k (l or l+1) of each state (in,integer(nst))
20 ! occ : occupancy of each state (inout,real(nst))
21 ! xctype : exchange-correlation type (in,integer(3))
22 ! xcgrad : 1 for GGA functional, 0 otherwise (in,integer)
23 ! nr : number of radial mesh points (in,integer)
24 ! r : radial mesh (in,real(nr))
25 ! eval : eigenvalue without rest-mass energy for each state (out,real(nst))
26 ! rho : charge density (out,real(nr))
27 ! vr : self-constistent potential (out,real(nr))
28 ! rwf : major and minor components of radial wavefunctions for each state
29 ! (out,real(nr,2,nst))
30 ! !DESCRIPTION:
31 ! Solves the Dirac-Kohn-Sham equations for an atom using the
32 ! exchange-correlation functional {\tt xctype} and returns the self-consistent
33 ! radial wavefunctions, eigenvalues, charge densities and potentials. Requires
34 ! the exchange-correlation interface routine {\tt xcifc}.
35 !
36 ! !REVISION HISTORY:
37 ! Created September 2002 (JKD)
38 ! Fixed s.c. convergence problem, October 2003 (JKD)
39 ! Added support for GGA functionals, June 2006 (JKD)
40 !
41 !EOP
42 !BOC
43 implicit none
44 ! arguments
45 real(8), intent(in) :: sol
46 logical, intent(in) :: ptnucl
47 real(8), intent(in) :: zn
48 integer, intent(in) :: nst
49 integer, intent(in) :: n(nst),l(nst),k(nst)
50 real(8), intent(inout) :: occ(nst)
51 integer, intent(in) :: xctype(3),xcgrad
52 integer, intent(in) :: nr
53 real(8), intent(in) :: r(nr)
54 real(8), intent(out) :: eval(nst)
55 real(8), intent(out) :: rho(nr),vr(nr)
56 real(8), intent(out) :: rwf(nr,2,nst)
57 ! local variables
58 integer, parameter :: maxscl=200
59 integer ir,ist,iscl
60 real(8), parameter :: fourpi=12.566370614359172954d0
61 ! potential convergence tolerance
62 real(8), parameter :: eps=1.d-6
63 real(8) dv,dvp,ze,beta,t1
64 ! allocatable arrays
65 real(8), allocatable :: vn(:),vh(:),ex(:),ec(:),vx(:),vc(:),vrp(:)
66 real(8), allocatable :: ri(:),wpr(:,:),fr1(:),fr2(:),gr1(:),gr2(:)
67 real(8), allocatable :: grho(:),g2rho(:),g3rho(:)
68 if (nst < 1) then
69  write(*,*)
70  write(*,'("Error(atom): nst < 1 : ",I8)') nst
71  write(*,*)
72  stop
73 end if
74 ! allocate local arrays
75 allocate(vn(nr),vh(nr),ex(nr),ec(nr),vx(nr),vc(nr),vrp(nr))
76 allocate(ri(nr),wpr(4,nr),fr1(nr),fr2(nr),gr1(nr),gr2(nr))
77 if (xcgrad == 1) then
78  allocate(grho(nr),g2rho(nr),g3rho(nr))
79 end if
80 ! find total electronic charge
81 ze=0.d0
82 do ist=1,nst
83  ze=ze+occ(ist)
84 end do
85 ! set up nuclear potential
86 call potnucl(ptnucl,nr,r,zn,vn)
87 do ir=1,nr
88  ri(ir)=1.d0/r(ir)
89 ! initialise the Kohn-Sham potential to the nuclear potential
90  vr(ir)=vn(ir)
91 end do
92 ! determine the weights for radial integration
93 call wsplintp(nr,r,wpr)
94 dvp=0.d0
95 vrp(1:nr)=0.d0
96 ! initialise mixing parameter
97 beta=0.5d0
98 ! initialise eigenvalues to relativistic values (minus the rest mass energy)
99 do ist=1,nst
100  t1=sqrt(dble(k(ist)**2)-(zn/sol)**2)
101  t1=(dble(n(ist)-abs(k(ist)))+t1)**2
102  t1=1.d0+((zn/sol)**2)/t1
103  eval(ist)=sol**2/sqrt(t1)-sol**2
104 end do
105 ! start of self-consistent loop
106 do iscl=1,maxscl
107 ! solve the Dirac equation for each state
108 !$OMP PARALLEL DO DEFAULT(SHARED) &
109 !$OMP SCHEDULE(DYNAMIC)
110  do ist=1,nst
111  call rdirac(sol,n(ist),l(ist),k(ist),nr,r,vr,eval(ist),rwf(:,1,ist), &
112  rwf(:,2,ist))
113  end do
114 !$OMP END PARALLEL DO
115 ! compute the charge density
116  do ir=1,nr
117  t1=sum(occ(1:nst)*(rwf(ir,1,1:nst)**2+rwf(ir,2,1:nst)**2))
118  fr1(ir)=t1
119  fr2(ir)=t1*ri(ir)
120  rho(ir)=(1.d0/fourpi)*t1*ri(ir)**2
121  end do
122  call splintwp(nr,wpr,fr1,gr1)
123  call splintwp(nr,wpr,fr2,gr2)
124 ! find the Hartree potential
125  t1=gr2(nr)
126  vh(1:nr)=gr1(1:nr)*ri(1:nr)+t1-gr2(1:nr)
127 ! normalise charge density and potential
128  t1=ze/gr1(nr)
129  rho(1:nr)=t1*rho(1:nr)
130  vh(1:nr)=t1*vh(1:nr)
131 ! compute the exchange-correlation energy and potential
132  if (xcgrad == 1) then
133 ! GGA functional
134 ! |∇ρ|
135  call fderiv(1,nr,r,rho,grho)
136 ! ∇²ρ
137  call fderiv(2,nr,r,rho,g2rho)
138  do ir=1,nr
139  g2rho(ir)=g2rho(ir)+2.d0*ri(ir)*grho(ir)
140  end do
141 ! approximate (∇ρ)⋅(∇|∇ρ|)
142  do ir=1,nr
143  g3rho(ir)=grho(ir)*g2rho(ir)
144  end do
145  call xcifc(xctype,nr,rho=rho,grho=grho,g2rho=g2rho,g3rho=g3rho,ex=ex,ec=ec,&
146  vx=vx,vc=vc)
147  else
148 ! LDA functional
149  call xcifc(xctype,nr,rho=rho,ex=ex,ec=ec,vx=vx,vc=vc)
150  end if
151 ! self-consistent potential
152  vr(1:nr)=vh(1:nr)+vx(1:nr)+vc(1:nr)
153 ! determine change in potential
154  t1=sum((vr(1:nr)-vrp(1:nr))**2)
155  dv=sqrt(t1)/dble(nr)
156  if (iscl > 2) then
157 ! reduce beta if change in potential is diverging
158  if (dv > dvp) beta=beta*0.8d0
159  beta=max(beta,0.01d0)
160  end if
161  dvp=dv
162  do ir=1,nr
163 ! mix old and new potentials
164  vr(ir)=(1.d0-beta)*vrp(ir)+beta*vr(ir)
165  vrp(ir)=vr(ir)
166 ! add nuclear potential
167  vr(ir)=vr(ir)+vn(ir)
168  end do
169 ! check for convergence
170  if ((iscl > 2).and.(dv < eps)) goto 10
171 ! end self-consistent loop
172 end do
173 write(*,*)
174 write(*,'("Warning(atom): maximum iterations exceeded")')
175 10 continue
176 deallocate(vn,vh,ex,ec,vx,vc,vrp)
177 deallocate(ri,wpr,fr1,fr2,gr1,gr2)
178 if (xcgrad == 1) deallocate(grho,g2rho,g3rho)
179 
180 contains
181 
182 pure subroutine splintwp(n,wp,f,g)
183 implicit none
184 ! arguments
185 integer, intent(in) :: n
186 real(8), intent(in) :: wp(*),f(n)
187 real(8), intent(out) :: g(n)
188 ! local variables
189 integer i,j
190 real(8) sm
191 g(1)=0.d0
192 sm=wp(5)*f(1)+wp(6)*f(2)+wp(7)*f(3)+wp(8)*f(4)
193 g(2)=sm
194 do i=2,n-2
195  j=i*4+1
196  sm=sm+wp(j)*f(i-1)+wp(j+1)*f(i)+wp(j+2)*f(i+1)+wp(j+3)*f(i+2)
197  g(i+1)=sm
198 end do
199 j=(n-1)*4+1
200 g(n)=sm+wp(j)*f(n-3)+wp(j+1)*f(n-2)+wp(j+2)*f(n-1)+wp(j+3)*f(n)
201 end subroutine
202 
203 end subroutine
204 !EOC
205 
subroutine potnucl(ptnucl, nr, r, zn, vn)
Definition: potnucl.f90:10
subroutine atom(sol, ptnucl, zn, nst, n, l, k, occ, xctype, xcgrad, nr, r, eval, rho, vr, rwf)
Definition: atom.f90:10
subroutine fderiv(m, n, x, f, g)
Definition: fderiv.f90:10
subroutine rdirac(sol, n, l, k, nr, r, vr, eval, g0, f0)
Definition: rdirac.f90:10
pure subroutine splintwp(n, wp, f, g)
Definition: atom.f90:183
subroutine wsplintp(n, x, wp)
Definition: wsplintp.f90:7
subroutine xcifc(xctype, n, c_tb09, tempa, rho, rhoup, rhodn, grho, gup, gdn, g2rho, g2up, g2dn, g3rho, g3up, g3dn, grho2, gup2, gdn2, gupdn, tau, tauup, taudn, ex, ec, vx, vc, vxup, vxdn, vcup, vcdn, dxdgr2, dxdgu2, dxdgd2, dxdgud, dcdgr2, dcdgu2, dcdgd2, dcdgud, dxdg2r, dxdg2u, dxdg2d, dcdg2r, dcdg2u, dcdg2d, wx, wxup, wxdn, wc, wcup, wcdn, dtdr, dtdru, dtdrd, dtdgr2, dtdgu2, dtdgd2, dtdg2r, dtdg2u, dtdg2d)
Definition: modxcifc.f90:20