2003-01-23 22:51:02 +08:00
|
|
|
!
|
2009-10-20 14:18:03 +08:00
|
|
|
! Copyright (C) 2001-2009 Quantum ESPRESSO group
|
2005-03-21 22:01:19 +08:00
|
|
|
! 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-01 01:55:33 +08:00
|
|
|
subroutine adddvepsi_us(becp2,ipol,kpoint)
|
|
|
|
! This subdoutine adds to dvpsi the terms which depend on the augmentation
|
2009-10-20 14:18:03 +08:00
|
|
|
! charge. It assumes that the variable dpqq, has been set and it is in
|
2009-10-20 20:50:28 +08:00
|
|
|
! the crystal basis.
|
2009-02-12 16:56:50 +08:00
|
|
|
! It calculates the last two terms of Eq.10 in JCP 21, 9934 (2004).
|
|
|
|
! P^+_c is applied in solve_e.
|
2004-06-01 01:55:33 +08:00
|
|
|
!
|
2003-01-23 22:51:02 +08:00
|
|
|
|
2008-08-24 01:55:06 +08:00
|
|
|
USE kinds, only : DP
|
|
|
|
USE spin_orb, ONLY : lspinorb
|
|
|
|
USE uspp, ONLY : nkb, vkb, qq, qq_so
|
|
|
|
USE wvfct, ONLY : npwx, npw, nbnd
|
|
|
|
USE cell_base, ONLY : at
|
2004-06-12 21:44:18 +08:00
|
|
|
USE ions_base, ONLY : nat, ityp, ntyp => nsp
|
2007-10-09 00:17:11 +08:00
|
|
|
USE noncollin_module, ONLY : noncolin, npol
|
2004-06-01 01:55:33 +08:00
|
|
|
USE uspp_param, only: nh
|
2009-09-21 21:38:34 +08:00
|
|
|
USE phus, ONLY : becp1, dpqq, dpqq_so
|
2009-10-20 20:50:28 +08:00
|
|
|
USE becmod, ONLY : bec_type
|
2009-02-02 18:52:58 +08:00
|
|
|
USE control_ph, ONLY: nbnd_occ
|
|
|
|
USE eqv, ONLY : dvpsi
|
2004-06-01 01:55:33 +08:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer, intent(in) :: ipol, kpoint
|
2009-10-20 20:50:28 +08:00
|
|
|
TYPE(bec_type), intent(in) :: becp2
|
2004-06-01 01:55:33 +08:00
|
|
|
|
2009-10-20 14:18:03 +08:00
|
|
|
complex(DP), allocatable :: ps(:), ps_nc(:,:)
|
2009-10-20 20:50:28 +08:00
|
|
|
integer:: ijkb0, nt, na, ih, jh, ikb, jkb, ibnd, ip, is, js, ijs
|
2004-06-01 01:55:33 +08:00
|
|
|
|
2007-10-09 00:17:11 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
allocate (ps_nc(nbnd,npol))
|
|
|
|
ELSE
|
|
|
|
allocate (ps(nbnd))
|
|
|
|
END IF
|
2004-06-01 01:55:33 +08:00
|
|
|
|
|
|
|
ijkb0 = 0
|
|
|
|
do nt = 1, ntyp
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp(na).eq.nt) then
|
|
|
|
do ih = 1, nh (nt)
|
|
|
|
ikb = ijkb0 + ih
|
2007-10-09 00:17:11 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
ps_nc = (0.d0,0.d0)
|
|
|
|
ELSE
|
|
|
|
ps = (0.d0,0.d0)
|
|
|
|
END IF
|
2004-06-01 01:55:33 +08:00
|
|
|
do jh = 1, nh (nt)
|
|
|
|
jkb = ijkb0 + jh
|
|
|
|
do ibnd=1, nbnd_occ(kpoint)
|
2007-10-09 00:17:11 +08:00
|
|
|
IF (noncolin) THEN
|
2009-10-20 20:50:28 +08:00
|
|
|
IF (lspinorb) THEN
|
|
|
|
ijs=0
|
|
|
|
do is=1,npol
|
|
|
|
do js=1,npol
|
|
|
|
ijs=ijs+1
|
|
|
|
ps_nc(ibnd,is)=ps_nc(ibnd,is) + &
|
|
|
|
qq_so(ih,jh,ijs,nt)* &
|
|
|
|
(0.d0,1.d0)*becp2%nc(jkb,js,ibnd) &
|
|
|
|
+ becp1(kpoint)%nc(jkb,js,ibnd)* &
|
|
|
|
dpqq_so(ih,jh,ijs,ipol,nt)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
ELSE
|
|
|
|
DO is=1,npol
|
|
|
|
ps_nc(ibnd,is)=ps_nc(ibnd,is)+ &
|
|
|
|
qq(ih,jh,nt)*becp2%nc(jkb,is,ibnd)*(0.d0,1.d0) &
|
|
|
|
+ dpqq(ih,jh,ipol,nt)* &
|
|
|
|
becp1(kpoint)%nc(jkb,is,ibnd)
|
|
|
|
END DO
|
|
|
|
END IF
|
2007-10-09 00:17:11 +08:00
|
|
|
ELSE
|
2009-10-20 20:50:28 +08:00
|
|
|
ps(ibnd) = ps(ibnd)+qq(ih,jh,nt)*becp2%k(jkb,ibnd) &
|
|
|
|
*(0.d0,1.d0) + &
|
|
|
|
dpqq(ih,jh,ipol,nt)* becp1(kpoint)%k(jkb,ibnd)
|
2007-10-09 00:17:11 +08:00
|
|
|
END IF
|
2004-06-01 01:55:33 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
do ibnd = 1, nbnd_occ (kpoint)
|
2007-10-09 00:17:11 +08:00
|
|
|
IF (noncolin) THEN
|
2009-08-03 17:19:02 +08:00
|
|
|
CALL zaxpy(npw,ps_nc(ibnd,1),vkb(1,ikb),1, &
|
2007-10-09 00:17:11 +08:00
|
|
|
dvpsi(1,ibnd),1)
|
2009-08-03 17:19:02 +08:00
|
|
|
CALL zaxpy(npw,ps_nc(ibnd,2),vkb(1,ikb),1, &
|
2007-10-09 00:17:11 +08:00
|
|
|
dvpsi(1+npwx,ibnd),1)
|
|
|
|
ELSE
|
2009-08-03 17:19:02 +08:00
|
|
|
CALL zaxpy(npw,ps(ibnd),vkb(1,ikb),1,dvpsi(1,ibnd),1)
|
2007-10-09 00:17:11 +08:00
|
|
|
END IF
|
2004-06-01 01:55:33 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
ijkb0=ijkb0+nh(nt)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
if (jkb.ne.nkb) call errore ('adddvepsi_us', 'unexpected error', 1)
|
|
|
|
|
2007-10-09 00:17:11 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
deallocate(ps_nc)
|
|
|
|
ELSE
|
|
|
|
deallocate(ps)
|
|
|
|
END IF
|
|
|
|
|
2004-06-01 01:55:33 +08:00
|
|
|
|
2007-10-09 00:17:11 +08:00
|
|
|
RETURN
|
|
|
|
END SUBROUTINE adddvepsi_us
|