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 .
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine compute_alphasum
|
2004-03-07 21:47:42 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine computes the alphasum term which is used to compute the
|
|
|
|
! change of the charge due to the displacement of the augmentation
|
|
|
|
! term. (See Eq. 29)
|
|
|
|
! It implements Eq.13 of the notes.
|
|
|
|
!
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-06-12 21:44:18 +08:00
|
|
|
|
|
|
|
USE ions_base, ONLY : nat, ityp, ntyp => nsp
|
2004-03-07 21:47:42 +08:00
|
|
|
use pwcom
|
2007-02-08 21:40:53 +08:00
|
|
|
USE noncollin_module, ONLY : noncolin, npol
|
2004-03-07 21:47:42 +08:00
|
|
|
USE kinds, only : DP
|
2004-06-01 01:55:33 +08:00
|
|
|
USE uspp_param, ONLY: nh, tvanp
|
2004-03-07 21:47:42 +08:00
|
|
|
use phcom
|
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
integer :: ik, ikk, ikq, ijkb0, ijh, ikb, jkb, ih, jh, na, nt, &
|
2007-02-08 21:40:53 +08:00
|
|
|
ipol, ibnd, is1, is2
|
2004-03-07 21:47:42 +08:00
|
|
|
! counter on k points
|
|
|
|
! counters on beta functions
|
|
|
|
! counters on beta functions
|
|
|
|
! counters for atoms
|
|
|
|
! counter on polarizations
|
|
|
|
! counter on bands
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: wgg1
|
2004-03-07 21:47:42 +08:00
|
|
|
! auxiliary weight
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (.not.okvan) return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2007-02-08 21:40:53 +08:00
|
|
|
alphasum = 0.d0
|
|
|
|
IF (noncolin) alphasum_nc=(0.d0,0.d0)
|
2004-03-07 21:47:42 +08:00
|
|
|
do ik = 1, nksq
|
|
|
|
if (lgamma) then
|
|
|
|
ikk = ik
|
|
|
|
ikq = ik
|
|
|
|
else
|
|
|
|
ikk = 2 * ik - 1
|
|
|
|
ikq = ikk + 1
|
|
|
|
endif
|
|
|
|
if (lsda) current_spin = isk (ikk)
|
|
|
|
ijkb0 = 0
|
|
|
|
do nt = 1, ntyp
|
|
|
|
if (tvanp (nt) ) then
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) == nt) then
|
|
|
|
ijh = 0
|
|
|
|
do ih = 1, nh (nt)
|
|
|
|
ikb = ijkb0 + ih
|
|
|
|
ijh = ijh + 1
|
|
|
|
do ibnd = 1, nbnd_occ (ikk)
|
|
|
|
wgg1 = wg (ibnd, ikk)
|
|
|
|
do ipol = 1, 3
|
2007-02-08 21:40:53 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
DO is1=1,npol
|
|
|
|
DO is2=1,npol
|
|
|
|
alphasum_nc(ijh,ipol,na,is1,is2) = &
|
|
|
|
alphasum_nc(ijh,ipol,na,is1,is2)+wgg1* &
|
|
|
|
(CONJG(alphap_nc(ikb,is1,ibnd,ipol,ik))*&
|
|
|
|
becp1_nc(ikb,is2,ibnd,ik) + &
|
|
|
|
CONJG(becp1_nc(ikb,is1,ibnd,ik))* &
|
|
|
|
alphap_nc(ikb,is2,ibnd,ipol,ik))
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
ELSE
|
|
|
|
alphasum(ijh,ipol,na,current_spin) = &
|
2004-03-07 21:47:42 +08:00
|
|
|
alphasum(ijh,ipol,na,current_spin) + 2.d0*wgg1*&
|
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
|
|
|
DBLE (CONJG(alphap (ikb,ibnd,ipol,ik) ) * &
|
2004-03-07 21:47:42 +08:00
|
|
|
becp1 (ikb,ibnd,ik) )
|
2007-02-08 21:40:53 +08:00
|
|
|
END IF
|
2004-03-07 21:47:42 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2007-02-08 21:40:53 +08:00
|
|
|
do jh = ih+1, nh (nt)
|
2004-03-07 21:47:42 +08:00
|
|
|
jkb = ijkb0 + jh
|
2007-02-08 21:40:53 +08:00
|
|
|
ijh = ijh + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
do ibnd = 1, nbnd
|
2007-02-08 21:40:53 +08:00
|
|
|
wgg1 = wg (ibnd, ikk)
|
|
|
|
do ipol = 1, 3
|
|
|
|
IF (noncolin) THEN
|
|
|
|
DO is1=1,npol
|
|
|
|
DO is2=1,npol
|
|
|
|
alphasum_nc(ijh,ipol,na,is1,is2) = &
|
|
|
|
alphasum_nc(ijh,ipol,na,is1,is2) &
|
|
|
|
+wgg1* &
|
|
|
|
(CONJG(alphap_nc(ikb,is1,ibnd,ipol,ik))* &
|
|
|
|
becp1_nc(jkb,is2,ibnd,ik)+ &
|
|
|
|
CONJG(becp1_nc(ikb,is1,ibnd,ik))* &
|
|
|
|
alphap_nc(jkb,is2,ibnd,ipol,ik) )
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
ELSE
|
|
|
|
alphasum(ijh,ipol,na,current_spin) = &
|
|
|
|
alphasum(ijh,ipol,na,current_spin) + &
|
|
|
|
2.d0 * wgg1 * &
|
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
|
|
|
DBLE (CONJG(alphap(ikb,ibnd,ipol,ik) )*&
|
2004-03-07 21:47:42 +08:00
|
|
|
becp1 (jkb,ibnd,ik) + &
|
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
|
|
|
CONJG( becp1 (ikb,ibnd,ik) ) * &
|
2004-03-07 21:47:42 +08:00
|
|
|
alphap (jkb,ibnd,ipol,ik) )
|
2007-02-08 21:40:53 +08:00
|
|
|
END IF
|
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
ijkb0 = ijkb0 + nh (nt)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
else
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) == nt) ijkb0 = ijkb0 + nh (nt)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
2007-02-08 21:40:53 +08:00
|
|
|
|
|
|
|
IF (noncolin.and.okvan) THEN
|
|
|
|
DO nt = 1, ntyp
|
|
|
|
IF ( tvanp(nt) ) THEN
|
|
|
|
DO na = 1, nat
|
|
|
|
IF (ityp(na)==nt) THEN
|
|
|
|
IF (so(nt)) THEN
|
|
|
|
CALL transform_alphasum_so(alphasum_nc,na)
|
|
|
|
ELSE
|
|
|
|
CALL transform_alphasum_nc(alphasum_nc,na)
|
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
END DO
|
|
|
|
END IF
|
|
|
|
END DO
|
|
|
|
END IF
|
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
! do na=1,nat
|
|
|
|
! nt=ityp(na)
|
|
|
|
! do ijh=1,nh(nt)*(nh(nt)+1)/2
|
|
|
|
! do ipol=1,3
|
|
|
|
! WRITE( stdout,'(3i5,f20.10)') na, ijh, ipol,
|
|
|
|
! + alphasum(ijh,ipol,na,1)
|
|
|
|
! enddo
|
|
|
|
! enddo
|
|
|
|
! enddo
|
|
|
|
! call stop_ph(.true.)
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine compute_alphasum
|