Further cleanup using becp_type.

git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@6057 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
dalcorso 2009-10-20 13:46:17 +00:00
parent 6d2912d456
commit 2905434654
1 changed files with 42 additions and 61 deletions

View File

@ -21,7 +21,7 @@ subroutine add_for_charges (ik, uact)
USE spin_orb, ONLY : lspinorb
USE uspp, ONLY : nkb, qq, qq_so, vkb
USE wvfct, ONLY : npwx, npw, nbnd, igk
USE becmod, ONLY: calbec
USE becmod, ONLY: calbec, bec_type, allocate_bec_type, deallocate_bec_type
USE noncollin_module, ONLY : noncolin, npol
USE uspp_param, only: nh
USE eqv, ONLY : dvpsi, dpsi
@ -42,7 +42,7 @@ subroutine add_for_charges (ik, uact)
!
integer :: na, nb, mu, nu, ikk, ikq, ig, igg, nt, ibnd, ijkb0, &
ikb, jkb, ih, jh, ipol, is
ikb, jkb, ih, jh, ipol, is, js, ijs
! counter on atoms
! counter on modes
! the point k
@ -65,22 +65,22 @@ subroutine add_for_charges (ik, uact)
! the scalar product
! the scalar product
! a mesh space for psi
complex(DP), allocatable :: bedp(:,:), aux1(:,:), alphapp(:,:,:)
complex(DP), allocatable :: bedp_nc(:,:,:), alphapp_nc(:,:,:,:)
TYPE(bec_type) :: bedp, alphapp(3)
complex(DP), allocatable :: aux1(:,:)
logical :: ok
! used to save time
allocate (aux ( npwx))
allocate (aux1( npwx*npol, nbnd))
CALL allocate_bec_type(nkb,nbnd,bedp)
DO ipol=1,3
CALL allocate_bec_type(nkb,nbnd,alphapp(ipol))
ENDDO
IF (noncolin) THEN
allocate (bedp_nc( nkb, npol, nbnd) )
allocate (alphapp_nc (nkb, npol, nbnd,3))
allocate (ps1_nc ( nkb, npol, nbnd))
allocate (ps2_nc ( nkb, npol, nbnd , 3))
ELSE
allocate (bedp( nkb, nbnd) )
allocate (alphapp (nkb,nbnd,3))
allocate (ps1 ( nkb , nbnd))
allocate (ps2 ( nkb , nbnd , 3))
ENDIF
@ -97,24 +97,24 @@ subroutine add_for_charges (ik, uact)
if (noncolin) then
ps1_nc = (0.d0, 0.d0)
ps2_nc = (0.d0, 0.d0)
alphapp_nc = (0.d0,0.d0)
bedp_nc = (0.d0,0.d0)
bedp%nc = (0.d0,0.d0)
DO ipol=1,3
alphapp(ipol)%nc = (0.d0,0.d0)
END DO
else
ps1 = (0.d0, 0.d0)
ps2 = (0.d0, 0.d0)
alphapp = (0.d0,0.d0)
bedp = (0.d0,0.d0)
bedp%k = (0.d0,0.d0)
DO ipol=1,3
alphapp(ipol)%k = (0.d0,0.d0)
END DO
endif
aux1 = (0.d0, 0.d0)
!
! first we calculate the products of the beta functions with dpsi
!
IF (noncolin) THEN
call calbec (npw, vkb, dpsi, bedp_nc)
ELSE
call calbec (npw, vkb, dpsi, bedp)
ENDIF
CALL calbec (npw, vkb, dpsi, bedp)
do ipol = 1, 3
aux1=(0.d0,0.d0)
do ibnd = 1, nbnd
@ -131,11 +131,7 @@ subroutine add_for_charges (ik, uact)
enddo
endif
enddo
if (noncolin) then
call calbec ( npw, vkb, aux1, alphapp_nc(:,:,:,ipol) )
else
call calbec ( npw, vkb, aux1, alphapp(:,:,ipol) )
endif
CALL calbec ( npw, vkb, aux1, alphapp(ipol) )
enddo
@ -155,55 +151,41 @@ subroutine add_for_charges (ik, uact)
do ibnd = 1, nbnd
if (noncolin) then
if (lspinorb) then
ps1_nc(ikb,1,ibnd)=ps1_nc(ikb,1,ibnd) + &
(qq_so (ih, jh, 1, nt) * &
alphapp_nc(jkb, 1, ibnd, ipol) + &
qq_so (ih, jh, 2, nt) * &
alphapp_nc(jkb, 2, ibnd, ipol) )* &
uact (mu + ipol)
ps1_nc(ikb,2,ibnd)=ps1_nc(ikb,2,ibnd) + &
(qq_so (ih, jh, 3, nt) * &
alphapp_nc(jkb, 1, ibnd, ipol) + &
qq_so (ih, jh, 4, nt) * &
alphapp_nc(jkb, 2, ibnd, ipol) )* &
uact (mu + ipol)
ps2_nc(ikb,1,ibnd,ipol)= &
ps2_nc(ikb,1,ibnd,ipol) + &
(qq_so (ih, jh, 1, nt) * &
bedp_nc (jkb, 1, ibnd) + &
qq_so (ih, jh, 2, nt) * &
bedp_nc (jkb, 2, ibnd) )*(0.d0,-1.d0) * &
uact (mu + ipol) * tpiba
ps2_nc(ikb,2,ibnd,ipol)= &
ps2_nc(ikb,2,ibnd,ipol) + &
(qq_so (ih, jh, 3, nt) * &
bedp_nc (jkb, 1, ibnd) + &
qq_so (ih, jh, 4, nt) * &
bedp_nc (jkb, 2, ibnd) )*(0.d0,-1.d0)* &
uact (mu + ipol) * tpiba
ijs=0
DO is=1,npol
DO js=1,npol
ijs=ijs+1
ps1_nc(ikb,is,ibnd)=ps1_nc(ikb,is,ibnd)+&
(qq_so (ih, jh, ijs, nt) * &
alphapp(ipol)%nc(jkb,js,ibnd))* &
uact (mu + ipol)
ps2_nc(ikb,is,ibnd,ipol)= &
ps2_nc(ikb,is,ibnd,ipol) + &
(qq_so (ih, jh, ijs, nt) * &
bedp%nc (jkb, js, ibnd))*(0.d0,-1.d0)* &
uact (mu + ipol) * tpiba
ENDDO
ENDDO
else
do is=1,npol
ps1_nc(ikb,is,ibnd)=ps1_nc(ikb,is,ibnd) + &
qq (ih, jh, nt) * &
alphapp_nc(jkb, is, ibnd, ipol) * &
alphapp(ipol)%nc(jkb, is, ibnd) * &
uact (mu + ipol)
ps2_nc(ikb,is,ibnd,ipol)= &
ps2_nc(ikb,is, ibnd, ipol) + &
qq (ih, jh, nt) * (0.d0, -1.d0) * &
bedp_nc (jkb, is, ibnd) * &
bedp%nc (jkb, is, ibnd) * &
uact (mu + ipol) * tpiba
end do
endif
else
ps1 (ikb, ibnd) = ps1 (ikb, ibnd) + &
qq (ih, jh, nt) * &
alphapp(jkb, ibnd, ipol) * &
qq (ih, jh, nt)*alphapp(ipol)%k(jkb, ibnd)* &
uact (mu + ipol)
ps2 (ikb, ibnd, ipol) = ps2 (ikb, ibnd, ipol) + &
qq (ih, jh, nt) * &
(0.d0, -1.d0) * &
bedp (jkb, ibnd) * &
uact (mu + ipol) * tpiba
qq (ih, jh, nt) * (0.d0, -1.d0) * &
bedp%k(jkb, ibnd) *uact (mu + ipol) * tpiba
endif
enddo
enddo
@ -270,17 +252,16 @@ subroutine add_for_charges (ik, uact)
deallocate (aux)
deallocate (aux1)
IF (noncolin) THEN
deallocate (bedp_nc)
deallocate (alphapp_nc)
deallocate (ps1_nc)
deallocate (ps2_nc)
ELSE
deallocate (bedp)
deallocate (alphapp)
deallocate (ps1)
deallocate (ps2)
END IF
CALL deallocate_bec_type(bedp)
DO ipol=1,3
CALL deallocate_bec_type(alphapp(ipol))
END DO
return
end subroutine add_for_charges