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 : ",I0)') 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=sum(occ(1:nst))
82 ! set up nuclear potential
83 call potnucl(ptnucl,nr,r,zn,vn)
84 ! initialise the Kohn-Sham potential to the nuclear potential
85 vr(1:nr)=vn(1:nr)
86 ! pre-calculate 1/r
87 ri(1:nr)=1.d0/r(1:nr)
88 ! determine the weights for radial integration
89 call wsplintp(nr,r,wpr)
90 dvp=0.d0
91 vrp(1:nr)=0.d0
92 ! initialise mixing parameter
93 beta=0.5d0
94 ! initialise eigenvalues to relativistic values (minus the rest mass energy)
95 do ist=1,nst
96  t1=sqrt(dble(k(ist)**2)-(zn/sol)**2)
97  t1=(dble(n(ist)-abs(k(ist)))+t1)**2
98  t1=1.d0+((zn/sol)**2)/t1
99  eval(ist)=sol**2/sqrt(t1)-sol**2
100 end do
101 ! start of self-consistent loop
102 do iscl=1,maxscl
103 ! solve the Dirac equation for each state
104 !$OMP PARALLEL DO DEFAULT(SHARED) &
105 !$OMP SCHEDULE(DYNAMIC)
106  do ist=1,nst
107  call rdirac(sol,n(ist),l(ist),k(ist),nr,r,vr,eval(ist),rwf(:,1,ist), &
108  rwf(:,2,ist))
109  end do
110 !$OMP END PARALLEL DO
111 ! compute the charge density
112  do ir=1,nr
113  t1=sum(occ(1:nst)*(rwf(ir,1,1:nst)**2+rwf(ir,2,1:nst)**2))
114  fr1(ir)=t1
115  fr2(ir)=t1*ri(ir)
116  rho(ir)=(t1/fourpi)*ri(ir)**2
117  end do
118  call splintwp(nr,wpr,fr1,gr1)
119  call splintwp(nr,wpr,fr2,gr2)
120 ! find the Hartree potential
121  t1=gr2(nr)
122  vh(1:nr)=gr1(1:nr)*ri(1:nr)+t1-gr2(1:nr)
123 ! normalise charge density and potential
124  t1=ze/gr1(nr)
125  rho(1:nr)=t1*rho(1:nr)
126  vh(1:nr)=t1*vh(1:nr)
127 ! compute the exchange-correlation energy and potential
128  if (xcgrad == 1) then
129 ! GGA functional
130 ! |∇ρ|
131  call fderiv(1,nr,r,rho,grho)
132 ! ∇²ρ
133  call fderiv(2,nr,r,rho,g2rho)
134  do ir=1,nr
135  g2rho(ir)=g2rho(ir)+2.d0*ri(ir)*grho(ir)
136  end do
137 ! approximate (∇ρ)⋅(∇|∇ρ|)
138  do ir=1,nr
139  g3rho(ir)=grho(ir)*g2rho(ir)
140  end do
141  call xcifc(xctype,nr,rho=rho,grho=grho,g2rho=g2rho,g3rho=g3rho,ex=ex,ec=ec,&
142  vx=vx,vc=vc)
143  else
144 ! LDA functional
145  call xcifc(xctype,nr,rho=rho,ex=ex,ec=ec,vx=vx,vc=vc)
146  end if
147 ! self-consistent potential
148  vr(1:nr)=vh(1:nr)+vx(1:nr)+vc(1:nr)
149 ! determine change in potential
150  t1=sum((vr(1:nr)-vrp(1:nr))**2)
151  dv=sqrt(t1)/dble(nr)
152  if (iscl > 2) then
153 ! reduce beta if change in potential is diverging
154  if (dv > dvp) beta=beta*0.8d0
155  beta=max(beta,0.01d0)
156  end if
157  dvp=dv
158  do ir=1,nr
159 ! mix old and new potentials
160  vr(ir)=(1.d0-beta)*vrp(ir)+beta*vr(ir)
161  vrp(ir)=vr(ir)
162 ! add nuclear potential
163  vr(ir)=vr(ir)+vn(ir)
164  end do
165 ! check for convergence
166  if ((iscl > 2).and.(dv < eps)) goto 10
167 ! end self-consistent loop
168 end do
169 write(*,*)
170 write(*,'("Warning(atom): maximum iterations exceeded")')
171 10 continue
172 deallocate(vn,vh,ex,ec,vx,vc,vrp)
173 deallocate(ri,wpr,fr1,fr2,gr1,gr2)
174 if (xcgrad == 1) deallocate(grho,g2rho,g3rho)
175 
176 contains
177 
178 pure subroutine splintwp(n,wp,f,g)
179 implicit none
180 ! arguments
181 integer, intent(in) :: n
182 real(8), intent(in) :: wp(*),f(n)
183 real(8), intent(out) :: g(n)
184 ! local variables
185 integer i,j
186 real(8) sm
187 g(1)=0.d0
188 sm=wp(5)*f(1)+wp(6)*f(2)+wp(7)*f(3)+wp(8)*f(4)
189 g(2)=sm
190 do i=2,n-2
191  j=i*4+1
192  sm=sm+wp(j)*f(i-1)+wp(j+1)*f(i)+wp(j+2)*f(i+1)+wp(j+3)*f(i+2)
193  g(i+1)=sm
194 end do
195 j=(n-1)*4+1
196 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)
197 end subroutine
198 
199 end subroutine
200 !EOC
201 
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:179
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