quantum-espresso/PWCOND/four.f90

240 lines
6.7 KiB
Fortran

!
! Copyright (C) 2003 A. Smogunov
! 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 .
!
subroutine four(w0, z0, dz, tblm, taunew, r, rab, betar)
!
! This routine computes the bidimensional fourier transform of the
! beta function. It has been implemented for s, p, d-orbitals.
!
! w0(z,g,m)=1/S * \int w(r) \exp{-ig r_\perp} dr_\perp
! where w(r) - beta function of the alpha's orbital.
!
! (see Gradshtein "Tables of integrals")
! For a fixed l it computes w0 for all m.
!
! The order of spherical harmonics used:
! s ;
! p_z, p_{-x}, p_{-y} ;
! d_{z^2-1}, d_{-xz}, d_{-yz}, d_{x^2-y^2}, d_{xy}
!
! input: tblm - array characterizing the orbital.
! taunew - coordinates and radius of the orbital.
! z0 - the initial z
! dz - the slab width
!
! output: w0(z, g, m), where
! z0< z <z0+dz
! g - 2D g-vector
!
USE kinds, ONLY: DP
USE constants, ONLY : tpi, fpi
USE radial_grids, only : ndmx
USE cell_base, ONLY : alat, tpiba
USE cond, ONLY : sarea, nz1, ngper, gper, ninsh, gnsh, ngpsh
implicit none
integer :: kz, ig, ign, igphi, &
indexr, iz, lb, ir, nmesh, nmeshs, tblm(4)
real(DP), parameter :: eps=1.d-8
complex(DP), parameter :: cim=(0.d0, 1.d0)
real(DP) :: gn, s1, s2, cs, sn, cs2, sn2, rz, dz1, zr, &
dr, z0, dz, bessj, taunew(4), r(ndmx), &
rab(ndmx), betar(ndmx)
real(DP), allocatable :: x1(:), x2(:), x3(:), x4(:)
real(DP), allocatable :: fx1(:), fx2(:), fx3(:), fx4(:), zsl(:)
complex(DP) :: w0(nz1, ngper, 5)
complex(DP), allocatable :: wadd(:,:)
allocate( x1(0:ndmx) )
allocate( x2(0:ndmx) )
allocate( x3(0:ndmx) )
allocate( x4(0:ndmx) )
allocate( fx1( nz1 ) )
allocate( fx2( nz1 ) )
allocate( fx3( nz1 ) )
allocate( fx4( nz1 ) )
allocate( zsl( nz1) )
allocate( wadd( nz1, ngper ) )
lb = tblm(3)
nmesh=indexr(taunew(4)*alat,ndmx,r)
dz1=dz/nz1
zsl(1)=(z0+dz1*0.5d0-taunew(3))*alat
do kz = 2, nz1
zsl(kz) = zsl(kz-1)+dz1*alat
enddo
ig=0
do ign=1, ngpsh
gn=gnsh(ign)
do kz=1, nz1
if (abs(zsl(kz))+eps.le.taunew(4)*alat) then
iz=indexr(zsl(kz),nmesh,r)
if ((nmesh-iz)/2*2.eq.nmesh-iz) then
nmeshs=nmesh
else
nmeshs=nmesh+1
endif
do ir=iz, nmeshs
rz=sqrt(r(ir)**2-zsl(kz)**2)
if (lb.eq.0) then
x1(ir)=betar(ir)*bessj(0,gn*rz)
elseif (lb.eq.1) then
x1(ir)=betar(ir)*bessj(1,gn*rz)/r(ir)*rz
x2(ir)=betar(ir)*bessj(0,gn*rz)/r(ir)
elseif (lb.eq.2) then
x1(ir)=betar(ir)*bessj(2,gn*rz)*rz**2/r(ir)**2
x2(ir)=betar(ir)*bessj(1,gn*rz)*rz/r(ir)**2
x3(ir)=betar(ir)*bessj(0,gn*rz)/r(ir)**2
x4(ir)=betar(ir)*bessj(0,gn*rz)
else
call errore ('four','ls not programmed ',1)
endif
enddo
call simpson(nmeshs-iz+1,x1(iz),rab(iz),fx1(kz))
if (iz.eq.1) then
dr=r(iz)
else
dr=r(iz)-r(iz-1)
endif
zr=r(iz)-abs(zsl(kz))
if (lb.eq.0) then
if (iz.eq.1) then
x1(iz-1)=betar(iz)-betar(iz)/dr*zr
else
x1(iz-1)=betar(iz)-(betar(iz)-betar(iz-1))/dr*zr
endif
fx1(kz)=fx1(kz)+(x1(iz-1)+x1(iz))*0.5d0*zr
else
fx1(kz)=fx1(kz)+x1(iz)*0.5d0*zr
call simpson(nmeshs-iz+1,x2(iz),rab(iz),fx2(kz))
endif
if (lb.eq.1) then
if(iz.eq.1) then
x2(iz-1)=0.d0
else
x2(iz-1)=(betar(iz)-(betar(iz)- &
betar(iz-1))/dr*zr)/abs(zsl(kz))
endif
fx2(kz)=fx2(kz)+(x2(iz-1)+x2(iz))*0.5d0*zr
endif
if (lb.eq.2) then
fx2(kz)=fx2(kz)+x2(iz)*0.5d0*zr
call simpson(nmeshs-iz+1,x3(iz),rab(iz),fx3(kz))
call simpson(nmeshs-iz+1,x4(iz),rab(iz),fx4(kz))
if(iz.eq.1) then
x3(iz-1)=0.d0
x4(iz-1)=0.d0
else
x3(iz-1)=(betar(iz)-(betar(iz)- &
betar(iz-1))/dr*zr)/abs(zsl(kz))**2
x4(iz-1)=betar(iz)-(betar(iz)- &
betar(iz-1))/dr*zr
endif
fx3(kz)=fx3(kz)+(x3(iz-1)+x3(iz))*0.5d0*zr
fx4(kz)=fx4(kz)+(x4(iz-1)+x4(iz))*0.5d0*zr
endif
else
fx1(kz)=0.d0
fx2(kz)=0.d0
fx3(kz)=0.d0
fx4(kz)=0.d0
endif
enddo
do igphi=1, ninsh(ign)
ig=ig+1
if (gn.gt.eps) then
cs=gper(1,ig)*tpiba/gn
sn=gper(2,ig)*tpiba/gn
else
cs=0.d0
sn=0.d0
endif
cs2=cs**2-sn**2
sn2=2*cs*sn
do kz=1, nz1
if (lb.eq.0) then
w0(kz,ig,1)=fx1(kz)
elseif (lb.eq.1) then
w0(kz,ig,2)=cs*fx1(kz)
w0(kz,ig,1)=fx2(kz)
w0(kz,ig,3)=sn*fx1(kz)
elseif (lb.eq.2) then
w0(kz,ig,5)=sn2*fx1(kz)
w0(kz,ig,2)=cs*fx2(kz)
w0(kz,ig,1)=fx3(kz)
w0(kz,ig,3)=sn*fx2(kz)
w0(kz,ig,4)=cs2*fx1(kz)
wadd(kz,ig)=fx4(kz)
endif
enddo
enddo
enddo
if (lb.eq.0) then
s1=tpi/sarea/sqrt(fpi)
elseif (lb.eq.1) then
s1=tpi/sarea*sqrt(3.d0/fpi)
elseif (lb.eq.2) then
s1=-tpi/2.d0/sarea*sqrt(15.d0/fpi)
s2=tpi/sarea*sqrt(5.d0/tpi/8.d0)
endif
do ig=1, ngper
do kz=1, nz1
if (lb.eq.0) then
w0(kz,ig,1)=s1*w0(kz,ig,1)
elseif (lb.eq.1) then
w0(kz,ig,2)=cim*s1*w0(kz,ig,2)
w0(kz,ig,1)=s1*zsl(kz)*w0(kz,ig,1)
w0(kz,ig,3)=cim*s1*w0(kz,ig,3)
elseif (lb.eq.2) then
w0(kz,ig,5)=s1*w0(kz,ig,5)
w0(kz,ig,2)=-2.d0*cim*s1*zsl(kz)*w0(kz,ig,2)
w0(kz,ig,1)=3.d0*zsl(kz)**2*s2*w0(kz,ig,1)-s2*wadd(kz,ig)
w0(kz,ig,3)=-2.d0*cim*s1*zsl(kz)*w0(kz,ig,3)
w0(kz,ig,4)=s1*w0(kz,ig,4)
endif
enddo
enddo
deallocate(x1)
deallocate(x2)
deallocate(x3)
deallocate(x4)
deallocate(fx1)
deallocate(fx2)
deallocate(fx3)
deallocate(fx4)
deallocate(zsl)
deallocate(wadd)
return
end subroutine four
function indexr(zz, ndim, r)
USE kinds, only : DP
implicit none
integer :: iz, ndim, indexr
real(DP) :: zz, r(ndim)
!
! abs(zz)<r(indexr)
!
iz = 1
do while(r(iz).le.abs(zz)+1.d-10)
iz=iz+1
enddo
indexr=iz
return
end function indexr