2004-05-14 23:33:08 +08:00
|
|
|
!
|
2005-03-21 22:33:57 +08:00
|
|
|
! Copyright (C) 2004 PWSCF group
|
|
|
|
! This file is distributed under the terms of the
|
|
|
|
! GNU General Public License. See the file `License'
|
|
|
|
! in the root directory of the present distribution,
|
|
|
|
! or http://www.gnu.org/copyleft/gpl.txt .
|
|
|
|
!
|
|
|
|
!
|
2004-05-14 23:33:08 +08:00
|
|
|
!---------------------------------------------------------------
|
2007-10-21 22:59:56 +08:00
|
|
|
subroutine ascheqps(nam,lam,jam,e0,mesh,ndm,grid,vpot,thresh,y,beta,ddd,&
|
|
|
|
qq,nbeta,nwfx,lls,jjs,ikk,nstop)
|
2004-12-10 23:33:00 +08:00
|
|
|
!---------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! numerical integration of the radial schroedinger equation for
|
|
|
|
! bound states in a local potential.
|
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! This routine works at fixed e
|
|
|
|
! It allows a nonlocal potential and an overlap
|
2004-12-10 23:33:00 +08:00
|
|
|
! operator. Therefore it can solve a general schroedinger
|
2007-10-21 22:59:56 +08:00
|
|
|
! equation necessary to solve a US pseudopotential
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-06-26 15:38:38 +08:00
|
|
|
use io_global, only : stdout
|
2004-12-10 23:33:00 +08:00
|
|
|
use kinds, only : DP
|
2007-10-21 22:59:56 +08:00
|
|
|
use radial_grids, only: radial_grid_type, series
|
2004-05-14 23:33:08 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
implicit none
|
2007-08-12 08:08:53 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
type(radial_grid_type), intent(in) :: grid
|
2004-12-10 23:33:00 +08:00
|
|
|
integer :: &
|
2007-10-21 22:59:56 +08:00
|
|
|
nam, &
|
|
|
|
lam, & ! l angular momentum
|
|
|
|
mesh,& ! size of radial mesh
|
|
|
|
ndm, & ! maximum radial mesh
|
2004-12-10 23:33:00 +08:00
|
|
|
nbeta,& ! number of beta function
|
|
|
|
nwfx, & ! maximum number of beta functions
|
2007-10-21 22:59:56 +08:00
|
|
|
itscf, & ! scf iteration
|
|
|
|
ndcr, & ! number of required nodes
|
|
|
|
n1, n2, & ! counters
|
|
|
|
ikl, & ! auxiliary
|
|
|
|
nstop, & ! used to check the behaviour of the routine
|
2004-12-10 23:33:00 +08:00
|
|
|
lls(nbeta),&! for each beta the angular momentum
|
|
|
|
ikk(nbeta) ! for each beta the point where it become zero
|
2004-05-14 23:33:08 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
real(DP) :: &
|
|
|
|
e0,e, & ! output eigenvalue
|
|
|
|
jam, & ! j angular momentum
|
2004-12-10 23:33:00 +08:00
|
|
|
vpot(mesh),&! the local potential
|
2007-10-21 22:59:56 +08:00
|
|
|
thresh, & ! precision of eigenvalue
|
|
|
|
y(mesh), & ! the output solution
|
2004-12-10 23:33:00 +08:00
|
|
|
jjs(nwfx), & ! the j angular momentum
|
2007-10-21 22:59:56 +08:00
|
|
|
beta(ndm,nwfx),& ! the beta functions
|
|
|
|
work(nbeta), & ! auxiliary space
|
2004-12-10 23:33:00 +08:00
|
|
|
ddd(nwfx,nwfx),qq(nwfx,nwfx) ! parameters for computing B_ij
|
|
|
|
!
|
|
|
|
! the local variables
|
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
real(DP) :: &
|
|
|
|
ddx12, & ! dx**2/12 used for Numerov integration
|
2004-12-10 23:33:00 +08:00
|
|
|
sqlhf, & ! the term for angular momentum in equation
|
2007-10-21 22:59:56 +08:00
|
|
|
xl1, x4l6, x6l12, x8l20,& ! used for starting series expansion
|
|
|
|
ze2, & ! possible coulomb term aroun the origin (set 0)
|
|
|
|
b(0:3), & ! coefficients of taylor expansion of potential
|
|
|
|
c1,c2,c3,c4,b0e, & ! auxiliary for expansion of wavefunction
|
|
|
|
rr1,rr2, & ! values of y in the first points+ auxiliary
|
|
|
|
eup,elw, & ! actual energy interval
|
|
|
|
ymx, & ! the maximum value of the function
|
|
|
|
rap, & ! the ratio between the number of nodes
|
|
|
|
fe,integ,dfe,de, &! auxiliary for numerov computation of e
|
|
|
|
eps, & ! the epsilon of the delta e
|
|
|
|
yln, xp, expn,& ! used to compute the tail of the solution
|
|
|
|
int_0_inf_dr ! integral function
|
2004-05-14 23:33:08 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
real(DP), allocatable :: &
|
|
|
|
fun(:), & ! integrand function
|
|
|
|
f(:), & ! the f function
|
|
|
|
el(:),c(:) ! auxiliary for inward integration
|
|
|
|
|
|
|
|
integer, parameter :: &
|
|
|
|
maxter=100 ! maximum number of iterations
|
2004-05-14 23:33:08 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
integer :: &
|
|
|
|
n, & ! counter on mesh points
|
|
|
|
iter,& ! counter on iteration
|
|
|
|
ik, & ! matching point
|
|
|
|
ns, & ! counter on beta functions
|
|
|
|
l1, & ! lam+1
|
|
|
|
nst, & ! used in the integration routine
|
|
|
|
npt, & ! number of points for energy intervals
|
|
|
|
ninter,& ! number of possible energy intervals
|
|
|
|
icountn,& ! counter on energy intervals
|
|
|
|
ierr, &
|
|
|
|
ncross,& ! actual number of nodes
|
|
|
|
nstart ! starting point for inward integration
|
2007-08-12 08:08:53 +08:00
|
|
|
|
2008-04-24 23:26:26 +08:00
|
|
|
logical, save :: first(0:10,0:10) = .true.
|
|
|
|
|
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
if (mesh.ne.grid%mesh) call errore('compute_solution','mesh dimension is not as expected',1)
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
|
|
|
! set up constants and allocate variables the
|
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
allocate(fun(mesh), stat=ierr)
|
|
|
|
allocate(f(mesh), stat=ierr)
|
|
|
|
allocate(el(mesh), stat=ierr)
|
|
|
|
allocate(c(mesh), stat=ierr)
|
|
|
|
|
2007-10-11 23:14:49 +08:00
|
|
|
nstop=0
|
2008-02-20 23:21:07 +08:00
|
|
|
nstart=0
|
2007-10-21 22:59:56 +08:00
|
|
|
e=e0
|
|
|
|
! write(6,*) 'entering ', nam,lam, e
|
|
|
|
eup=0.3_DP*e
|
|
|
|
elw=1.3_dp*e
|
|
|
|
ndcr=nam-lam-1
|
|
|
|
! write(6,*) 'entering ascheqps', vpot(mesh-20)*grid%r(mesh-20)
|
|
|
|
|
|
|
|
ddx12=grid%dx*grid%dx/12.0_dp
|
|
|
|
l1=lam+1
|
|
|
|
nst=l1*2
|
|
|
|
sqlhf=(DBLE(lam)+0.5_dp)**2
|
|
|
|
xl1=lam+1
|
|
|
|
x4l6=4*lam+6
|
|
|
|
x6l12=6*lam+12
|
|
|
|
x8l20=8*lam+20
|
|
|
|
!
|
|
|
|
! series developement of the potential near the origin
|
|
|
|
!
|
|
|
|
do n=1,4
|
|
|
|
y(n)=vpot(n)
|
2004-12-10 23:33:00 +08:00
|
|
|
enddo
|
2007-10-21 22:59:56 +08:00
|
|
|
call series(y,grid%r,grid%r2,b)
|
|
|
|
|
|
|
|
! write(stdout,*) 'eneter lam,eup,elw,e',lam,nbeta,eup,elw,e
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! set up the f-function and determine the position of its last
|
|
|
|
! change of sign
|
|
|
|
! f < 0 (approximatively) means classically allowed region
|
|
|
|
! f > 0 " " " forbidden "
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
do iter=1,maxter
|
|
|
|
! write(6,*) 'starting iter', iter, elw, e, eup
|
|
|
|
ik=1
|
|
|
|
f(1)=ddx12*(grid%r2(1)*(vpot(1)-e)+sqlhf)
|
|
|
|
do n=2,mesh
|
|
|
|
f(n)=ddx12*(grid%r2(n)*(vpot(n)-e)+sqlhf)
|
|
|
|
if( f(n).ne.sign(f(n),f(n-1)).and.n.lt.mesh-5 ) ik=n
|
|
|
|
enddo
|
|
|
|
! if (ik.eq.0.and.nbeta.eq.0) ik=mesh*3/4
|
|
|
|
if (ik.eq.1.or.grid%r(ik)>4.0_DP) ik=mesh*3/4
|
|
|
|
|
|
|
|
if(ik.ge.mesh-2) then
|
2004-12-10 23:33:00 +08:00
|
|
|
do n=1,mesh
|
2007-10-21 22:59:56 +08:00
|
|
|
write(stdout,*) grid%r(n), vpot(n), f(n)
|
2004-12-10 23:33:00 +08:00
|
|
|
enddo
|
2007-10-21 22:59:56 +08:00
|
|
|
call errore('compute_solution', 'No point found for matching',1)
|
2004-12-10 23:33:00 +08:00
|
|
|
endif
|
2007-10-21 22:59:56 +08:00
|
|
|
!
|
|
|
|
! determine if ik is sufficiently large
|
|
|
|
!
|
|
|
|
do ns=1,nbeta
|
|
|
|
if (lls(ns).eq.lam .and. jjs(ns).eq.jam .and. ikk(ns).gt.ik) ik=ikk(ns)+3
|
2004-12-10 23:33:00 +08:00
|
|
|
enddo
|
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! if everything is ok continue the integration and define f
|
|
|
|
!
|
|
|
|
do n=1,mesh
|
|
|
|
f(n)=1.0_dp-f(n)
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! determination of the wave-function in the first two points by
|
|
|
|
! series developement
|
|
|
|
!
|
2007-12-12 06:46:03 +08:00
|
|
|
! no coulomb divergence in the origin for a pseudopotential
|
|
|
|
ze2=0.d0
|
2007-10-21 22:59:56 +08:00
|
|
|
b0e=b(0)-e
|
|
|
|
c1=0.5_dp*ze2/xl1
|
|
|
|
c2=(c1*ze2+b0e)/x4l6
|
|
|
|
c3=(c2*ze2+c1*b0e+b(1))/x6l12
|
|
|
|
c4=(c3*ze2+c2*b0e+c1*b(1)+b(2))/x8l20
|
|
|
|
!
|
|
|
|
rr1=(1.0_dp+grid%r(1)*(c1+grid%r(1)*(c2+grid%r(1)*(c3+grid%r(1)*c4))))*grid%r(1)**l1
|
|
|
|
rr2=(1.0_dp+grid%r(2)*(c1+grid%r(2)*(c2+grid%r(2)*(c3+grid%r(2)*c4))))*grid%r(2)**l1
|
|
|
|
y(1)=rr1/grid%sqr(1)
|
|
|
|
y(2)=rr2/grid%sqr(2)
|
|
|
|
!
|
|
|
|
! outward integration before ik
|
|
|
|
!
|
|
|
|
call integrate_outward (lam,jam,e,mesh,ndm,grid,f,b,y,beta,ddd,qq,&
|
|
|
|
nbeta,nwfx,lls,jjs,ikk,ik)
|
2004-05-14 23:33:08 +08:00
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
ncross=0
|
|
|
|
ymx=0.0_dp
|
|
|
|
do n=2,ik-1
|
|
|
|
if ( y(n) .ne. sign(y(n),y(n+1)) ) ncross=ncross+1
|
|
|
|
ymx=max(ymx,abs(y(n+1)))
|
|
|
|
end do
|
|
|
|
!
|
|
|
|
! If at this point the number of nodes is wrong it means that something
|
|
|
|
! is probably wrong in the calling routines. A ghost might be present
|
|
|
|
! in the pseudopotential. With a nonlocal pseudopotential there is no
|
|
|
|
! node theorem so strictly speaking the following instructions are
|
|
|
|
! wrong but sometimes they help so we keep them here.
|
|
|
|
!
|
2008-04-24 23:26:26 +08:00
|
|
|
if(ndcr /= ncross .and. first(nam,lam)) then
|
|
|
|
write(stdout,'(/,7x,5(a,i3))') 'WARNING! Expected number of nodes: ',ndcr, '= ',nam,'-',lam,&
|
|
|
|
'- 1, number of nodes found:', ncross,'.'
|
|
|
|
write(stdout,'(7x,a,/,7xa,/)') 'Setting wfc to zero for this iteration.',&
|
|
|
|
'(This warning will only be printed once per wavefunction)'
|
|
|
|
first(nam,lam) = .false.
|
|
|
|
endif
|
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
if(ndcr < ncross) then
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! too many crossings. e is an upper bound to the true eigen-
|
|
|
|
! value. increase abs(e)
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
|
|
|
|
eup=e
|
|
|
|
e=0.9_dp*elw+0.1_dp*eup
|
|
|
|
! write(6,*) 'too many crossing', ncross, ndcr
|
|
|
|
! call errore('aschqps','wrong number of nodes. Probably a Ghost?',1)
|
|
|
|
y=0.0_DP
|
|
|
|
ymx=0.0_dp
|
|
|
|
go to 300
|
|
|
|
else if (ndcr > ncross) then
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! too few crossings. e is a lower bound to the true eigen-
|
|
|
|
! value. decrease abs(e)
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
elw=e
|
|
|
|
e=0.9_dp*eup+0.1_dp*elw
|
|
|
|
! write(6,*) 'too few crossing', ncross, ndcr
|
|
|
|
! call errore('aschqps','wrong number of nodes. Probably a Ghost?',1)
|
|
|
|
y=0.0_DP
|
|
|
|
ymx=0.0_dp
|
|
|
|
go to 300
|
|
|
|
end if
|
|
|
|
!
|
|
|
|
! inward integration up to ik
|
|
|
|
!
|
|
|
|
call integrate_inward(e,mesh,ndm,grid,f,y,c,el,ik,nstart)
|
|
|
|
!
|
|
|
|
! if necessary, improve the trial eigenvalue by the cooley's
|
|
|
|
! procedure. jw cooley math of comp 15,363(1961)
|
|
|
|
!
|
|
|
|
fe=(12.0_dp-10.0_dp*f(ik))*y(ik)-f(ik-1)*y(ik-1)-f(ik+1)*y(ik+1)
|
|
|
|
!
|
|
|
|
! audjust the normalization if needed
|
|
|
|
!
|
|
|
|
if(ymx.ge.1.0e10_dp) y=y/ymx
|
|
|
|
!
|
|
|
|
! calculate the normalization
|
|
|
|
!
|
|
|
|
do n1=1,nbeta
|
|
|
|
if (lam.eq.lls(n1).and.abs(jam-jjs(n1)).lt.1.e-7_dp) then
|
|
|
|
ikl=ikk(n1)
|
|
|
|
do n=1,ikl
|
|
|
|
fun(n)=beta(n,n1)*y(n)*grid%sqr(n)
|
|
|
|
enddo
|
|
|
|
work(n1)=int_0_inf_dr(fun,grid,ikl,nst)
|
2004-12-10 23:33:00 +08:00
|
|
|
else
|
2007-10-21 22:59:56 +08:00
|
|
|
work(n1)=0.0_dp
|
2004-12-10 23:33:00 +08:00
|
|
|
endif
|
|
|
|
enddo
|
2007-10-21 22:59:56 +08:00
|
|
|
do n=1,nstart
|
|
|
|
fun(n)= y(n)*y(n)*grid%r(n)
|
|
|
|
enddo
|
|
|
|
integ=int_0_inf_dr(fun,grid,nstart,nst)
|
|
|
|
do n1=1,nbeta
|
|
|
|
do n2=1,nbeta
|
|
|
|
integ = integ + qq(n1,n2)*work(n1)*work(n2)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
dfe=-y(ik)*f(ik)/grid%dx/integ
|
|
|
|
de=-fe*dfe
|
|
|
|
eps=abs(de/e)
|
|
|
|
! write(6,'(i5, 3f20.12)') iter, e, de
|
|
|
|
if(abs(de).lt.thresh) go to 600
|
|
|
|
if(eps.gt.0.25_dp) de=0.25_dp*de/eps
|
|
|
|
if(de.gt.0.0_dp) elw=e
|
|
|
|
if(de.lt.0.0_dp) eup=e
|
|
|
|
e=e+de
|
|
|
|
if(e.gt.eup) e=0.9_dp*eup+0.1_dp*elw
|
|
|
|
if(e.lt.elw) e=0.9_dp*elw+0.1_dp*eup
|
|
|
|
300 continue
|
|
|
|
enddo
|
|
|
|
nstop=1
|
2008-02-20 23:21:07 +08:00
|
|
|
if (nstart==0) goto 900
|
2007-10-21 22:59:56 +08:00
|
|
|
600 continue
|
|
|
|
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! exponential tail of the solution if it was not computed
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
if (nstart.lt.mesh) then
|
|
|
|
do n=nstart,mesh-1
|
|
|
|
if (y(n) == 0.0_dp) then
|
|
|
|
y(n+1)=0.0_dp
|
|
|
|
else
|
|
|
|
yln=log(abs(y(n)))
|
|
|
|
xp=-sqrt(12.0_dp*abs(1.0_dp-f(n)))
|
|
|
|
expn=yln+xp
|
|
|
|
if (expn.lt.-80.0_dp) then
|
|
|
|
y(n+1)=0.0_dp
|
|
|
|
else
|
|
|
|
y(n+1)=sign(exp(expn),y(n))
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo
|
2004-12-10 23:33:00 +08:00
|
|
|
endif
|
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
! normalize the eigenfunction as if they were norm conserving.
|
|
|
|
! If this is a US PP the correct normalization is done outside this
|
|
|
|
! routine.
|
2004-12-10 23:33:00 +08:00
|
|
|
!
|
2007-10-21 22:59:56 +08:00
|
|
|
do n=1,mesh
|
|
|
|
el(n)=grid%r(n)*y(n)*y(n)
|
2004-12-10 23:33:00 +08:00
|
|
|
enddo
|
2007-10-21 22:59:56 +08:00
|
|
|
integ=int_0_inf_dr(el,grid,mesh,nst)
|
|
|
|
if (integ>0.0_DP) then
|
|
|
|
integ=sqrt(integ)
|
|
|
|
do n=1,mesh
|
|
|
|
y(n)=grid%sqr(n)*y(n)/integ
|
|
|
|
enddo
|
|
|
|
e0=e
|
2004-12-10 23:33:00 +08:00
|
|
|
else
|
2007-10-21 22:59:56 +08:00
|
|
|
nstop=1
|
2004-12-10 23:33:00 +08:00
|
|
|
endif
|
2007-10-21 22:59:56 +08:00
|
|
|
|
2007-12-19 22:07:33 +08:00
|
|
|
900 continue
|
|
|
|
|
2007-10-21 22:59:56 +08:00
|
|
|
deallocate(el)
|
|
|
|
deallocate(f )
|
|
|
|
deallocate(c )
|
|
|
|
deallocate(fun )
|
2004-12-10 23:33:00 +08:00
|
|
|
return
|
2007-10-21 22:59:56 +08:00
|
|
|
|
2004-12-10 23:33:00 +08:00
|
|
|
end subroutine ascheqps
|