2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2001 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-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-06-12 21:44:18 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine gen_us_dj (ik, dvkb)
|
2003-01-20 05:58:50 +08:00
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! Calculates the beta function pseudopotentials with
|
|
|
|
! the derivative of the Bessel functions
|
|
|
|
!
|
2004-06-12 21:44:18 +08:00
|
|
|
USE kinds, ONLY : DP
|
|
|
|
USE constants, ONLY : tpi
|
|
|
|
USE ions_base, ONLY : nat, ntyp => nsp, ityp, tau
|
|
|
|
USE cell_base, ONLY : tpiba
|
|
|
|
USE klist, ONLY : xk
|
|
|
|
USE gvect, ONLY : ig1, ig2, ig3, eigts1, eigts2, eigts3, g
|
|
|
|
USE wvfct, ONLY : npw, npwx, igk
|
|
|
|
USE uspp, ONLY : nkb, indv, nhtol, nhtolm
|
2007-07-13 17:43:45 +08:00
|
|
|
USE us, ONLY : nqx, tab, tab_d2y, dq, spline_ps
|
2006-09-15 17:06:15 +08:00
|
|
|
USE splinelib
|
2007-10-06 16:23:39 +08:00
|
|
|
USE uspp_param, ONLY : upf, lmaxkb, nbetam, nh
|
2004-06-12 21:44:18 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
implicit none
|
|
|
|
!
|
|
|
|
integer :: ik
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: dvkb (npwx, nkb)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! local variables
|
|
|
|
!
|
2004-01-22 20:44:53 +08:00
|
|
|
integer :: ikb, nb, ih, ig, i0, i1, i2, i3 , nt
|
2003-01-20 05:58:50 +08:00
|
|
|
! counter on beta functions
|
|
|
|
! counter on beta functions
|
|
|
|
! counter on beta functions
|
|
|
|
! counter on G vectors
|
|
|
|
! index of the first nonzero point in the r
|
|
|
|
! counter on atomic type
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: arg, px, ux, vx, wx
|
2003-01-20 05:58:50 +08:00
|
|
|
! argument of the atomic phase factor
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: phase, pref
|
2003-01-20 05:58:50 +08:00
|
|
|
! atomic phase factor
|
|
|
|
! prefactor
|
|
|
|
|
2006-12-14 16:53:47 +08:00
|
|
|
integer :: na, l, iig, lm
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP), allocatable :: djl (:,:,:), ylm (:,:), q (:), gk (:,:)
|
2006-12-14 16:53:47 +08:00
|
|
|
real(DP) :: qt, eps
|
2006-12-03 01:58:50 +08:00
|
|
|
parameter (eps = 1.0d-8)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: sk (:)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2007-01-25 19:42:52 +08:00
|
|
|
integer :: iq
|
2006-09-15 17:06:15 +08:00
|
|
|
real(DP), allocatable :: xdata(:)
|
|
|
|
|
2006-11-11 01:23:47 +08:00
|
|
|
if (nkb.eq.0) return
|
|
|
|
|
2004-01-22 20:44:53 +08:00
|
|
|
call start_clock('stres_us31')
|
|
|
|
|
2006-03-03 22:19:17 +08:00
|
|
|
allocate (djl( npw , nbetam , ntyp))
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (ylm( npw ,(lmaxkb + 1) **2))
|
|
|
|
allocate (gk( 3, npw))
|
|
|
|
allocate (q( npw))
|
|
|
|
do ig = 1, npw
|
2003-01-20 05:58:50 +08:00
|
|
|
gk (1,ig) = xk (1, ik) + g(1, igk(ig) )
|
|
|
|
gk (2,ig) = xk (2, ik) + g(2, igk(ig) )
|
|
|
|
gk (3,ig) = xk (3, ik) + g(3, igk(ig) )
|
|
|
|
q (ig) = gk(1, ig)**2 + gk(2, ig)**2 + gk(3, ig)**2
|
|
|
|
enddo
|
|
|
|
|
2004-01-22 20:44:53 +08:00
|
|
|
call stop_clock('stres_us31')
|
|
|
|
call start_clock('stres_us32')
|
2003-02-08 00:04:36 +08:00
|
|
|
call ylmr2 ((lmaxkb+1)**2, npw, gk, q, ylm)
|
2004-01-22 20:44:53 +08:00
|
|
|
call stop_clock('stres_us32')
|
|
|
|
call start_clock('stres_us33')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2007-01-25 19:42:52 +08:00
|
|
|
|
|
|
|
if (spline_ps) then
|
2007-07-13 17:43:45 +08:00
|
|
|
allocate(xdata(nqx))
|
|
|
|
do iq = 1, nqx
|
2007-01-25 19:42:52 +08:00
|
|
|
xdata(iq) = (iq - 1) * dq
|
|
|
|
enddo
|
|
|
|
endif
|
2006-09-15 17:06:15 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
do nt = 1, ntyp
|
2007-10-06 16:23:39 +08:00
|
|
|
do nb = 1, upf(nt)%nbeta
|
2003-02-08 00:04:36 +08:00
|
|
|
do ig = 1, npw
|
|
|
|
qt = sqrt(q (ig)) * tpiba
|
2007-01-25 19:42:52 +08:00
|
|
|
if (spline_ps) then
|
|
|
|
djl(ig,nb,nt) = splint_deriv(xdata, tab(:,nb,nt), &
|
|
|
|
tab_d2y(:,nb,nt), qt)
|
|
|
|
else
|
|
|
|
px = qt / dq - int (qt / dq)
|
|
|
|
ux = 1.d0 - px
|
|
|
|
vx = 2.d0 - px
|
|
|
|
wx = 3.d0 - px
|
|
|
|
i0 = qt / dq + 1
|
|
|
|
i1 = i0 + 1
|
|
|
|
i2 = i0 + 2
|
|
|
|
i3 = i0 + 3
|
|
|
|
djl(ig,nb,nt) = ( tab (i0, nb, nt) * (-vx*wx-ux*wx-ux*vx)/6.d0 + &
|
|
|
|
tab (i1, nb, nt) * (+vx*wx-px*wx-px*vx)/2.d0 - &
|
|
|
|
tab (i2, nb, nt) * (+ux*wx-px*wx-px*ux)/2.d0 + &
|
|
|
|
tab (i3, nb, nt) * (+ux*vx-px*vx-px*ux)/6.d0 )/dq
|
|
|
|
endif
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2004-01-22 20:44:53 +08:00
|
|
|
call stop_clock('stres_us33')
|
|
|
|
call start_clock('stres_us34')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (q)
|
|
|
|
deallocate (gk)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (sk( npw))
|
|
|
|
ikb = 0
|
|
|
|
do nt = 1, ntyp
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) .eq.nt) then
|
2003-01-20 05:58:50 +08:00
|
|
|
arg = (xk (1, ik) * tau(1,na) + &
|
|
|
|
xk (2, ik) * tau(2,na) + &
|
|
|
|
xk (3, ik) * tau(3,na) ) * tpi
|
General cleanup of intrinsic functions:
conversion to real => DBLE
(including real part of a complex number)
conversion to complex => CMPLX
complex conjugate => CONJG
imaginary part => AIMAG
All functions are uppercase.
CMPLX is preprocessed by f_defs.h and performs an explicit cast:
#define CMPLX(a,b) cmplx(a,b,kind=DP)
This implies that 1) f_defs.h must be included whenever a CMPLX is present,
2) CMPLX should stay in a single line, 3) DP must be defined.
All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx
removed - please do not reintroduce any of them.
Tested only with ifc7 and g95 - beware unintended side effects
Maybe not the best solution (explicit casts everywhere would be better)
but it can be easily changed with a script if the need arises.
The following code might be used to test for possible trouble:
program test_intrinsic
implicit none
integer, parameter :: dp = selected_real_kind(14,200)
real (kind=dp) :: a = 0.123456789012345_dp
real (kind=dp) :: b = 0.987654321098765_dp
complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp)
print *, ' A = ', a
print *, ' DBLE(A)= ', DBLE(a)
print *, ' C = ', c
print *, 'CONJG(C)= ', CONJG(c)
print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c)
print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp)
end program test_intrinsic
Note that CMPLX and REAL without a cast yield single precision numbers on
ifc7 and g95 !!!
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
|
|
|
phase = CMPLX (cos (arg), - sin (arg) )
|
2003-02-08 00:04:36 +08:00
|
|
|
do ig = 1, npw
|
|
|
|
iig = igk (ig)
|
2003-01-20 05:58:50 +08:00
|
|
|
sk (ig) = eigts1 (ig1 (iig), na) * &
|
|
|
|
eigts2 (ig2 (iig), na) * &
|
|
|
|
eigts3 (ig3 (iig), na) * phase
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do ih = 1, nh (nt)
|
|
|
|
nb = indv (ih, nt)
|
|
|
|
l = nhtol (ih, nt)
|
2004-04-28 18:25:36 +08:00
|
|
|
lm= nhtolm(ih, nt)
|
2003-02-08 00:04:36 +08:00
|
|
|
ikb = ikb + 1
|
2004-05-06 21:06:16 +08:00
|
|
|
pref = (0.d0, -1.d0) **l
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do ig = 1, npw
|
2003-01-20 05:58:50 +08:00
|
|
|
dvkb (ig, ikb) = djl (ig, nb, nt) * sk (ig) * ylm (ig, lm) &
|
|
|
|
* pref
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
2004-01-22 20:44:53 +08:00
|
|
|
call stop_clock('stres_us34')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-21 22:57:00 +08:00
|
|
|
if (ikb.ne.nkb) call errore ('gen_us_dj', 'unexpected error', 1)
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (sk)
|
|
|
|
deallocate (ylm)
|
|
|
|
deallocate (djl)
|
2009-01-30 22:41:13 +08:00
|
|
|
if (spline_ps) deallocate(xdata)
|
2003-02-08 00:04:36 +08:00
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine gen_us_dj
|
|
|
|
|