Implementation of LDA+U gamma_only calculations.

git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@588 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
fabris 2004-02-12 12:32:06 +00:00
parent 6bb7f06ec2
commit 73d33ae8ad
5 changed files with 84 additions and 37 deletions

View File

@ -81,6 +81,7 @@ SUBROUTINE forces
!
! ... The Hubbard contribution
!
IF ( lda_plus_u .AND. gamma_only) STOP 'gamma-only calculation of forces for LDA+U is not implemented'
IF ( lda_plus_u ) CALL force_hub( forceh )
!
! ... The ionic contribution is computed here

View File

@ -25,9 +25,10 @@ subroutine new_ns
Hubbard_alpha, swfcatom, eth, d1, d2, d3
USE lsda_mod, ONLY: lsda, current_spin, nspin, isk
USE symme, ONLY: nsym, irt
USE wvfct, ONLY: nbnd, npw, npwx, igk, wg
USE wvfct, ONLY: nbnd, npw, npwx, igk, wg, gamma_only
USE wavefunctions_module, ONLY : evc
USE us, ONLY: newpseudo
USE gvect, ONLY : gstart
use io_files
#ifdef __PARA
use para
@ -45,6 +46,7 @@ subroutine new_ns
real(kind=DP) :: t0, scnds
! cpu time spent
real(kind=DP), external :: DDOT
complex(kind=DP) :: ZDOTC
complex(kind=DP) , allocatable :: proj(:,:)
@ -82,19 +84,28 @@ subroutine new_ns
if (lsda) current_spin = isk(ik)
if (nks.gt.1) read (iunigk) npw, igk
if (nks.gt.1) call davcio (evc, nwordwfc, iunwfc, ik, - 1)
call davcio (swfcatom, nwordatwfc, iunat, ik, - 1)
!
! make the projection
!
do ibnd = 1, nbnd
do i = 1, natomwfc
proj (i, ibnd) = ZDOTC (npw, swfcatom (1, i), 1, evc (1, ibnd), 1)
do ibnd = 1, nbnd
do i = 1, natomwfc
if ( gamma_only ) then
proj (i, ibnd) = 2.d0 * &
DDOT(2*npw, swfcatom (1, i), 1, evc (1, ibnd), 1)
if (gstart.eq.2) proj (i, ibnd) = proj (i, ibnd) - &
swfcatom (1, i) * evc (1, ibnd)
else
proj (i, ibnd) = ZDOTC (npw, swfcatom (1, i), 1, evc (1, ibnd), 1)
endif
enddo
enddo
enddo
#ifdef __PARA
call reduce (2 * natomwfc * nbnd, proj)
call reduce (2 * natomwfc * nbnd, proj)
#endif
!
! compute the occupation numbers (the quantities n(m1,m2)) of the
! atomic orbitals
@ -106,15 +117,14 @@ subroutine new_ns
do m2 = m1, 2 * Hubbard_l(nt) + 1
do ibnd = 1, nbnd
nr(m1,m2,current_spin,na) = nr(m1,m2,current_spin,na) + &
wg(ibnd,ik) * DREAL( proj(offset(na)+m2,ibnd) * &
conjg(proj(offset(na)+m1,ibnd)) )
wg(ibnd,ik) * DREAL( proj(offset(na)+m2,ibnd) * &
conjg(proj(offset(na)+m1,ibnd)) )
enddo
enddo
enddo
endif
enddo
! on k-points
enddo
! on k-points
enddo
#ifdef __PARA
@ -196,8 +206,9 @@ subroutine new_ns
enddo
enddo
enddo
endif
endif
enddo
!
! Now the contribution to the total energy is computed. The corrections
! needed to obtain a variational expression are already included
@ -216,6 +227,7 @@ subroutine new_ns
enddo
endif
enddo
deallocate ( offset, proj, nr )
if (nspin.eq.1) eth = 2.d0 * eth

View File

@ -22,13 +22,17 @@ subroutine orthoatwfc
USE basis, ONLY: nat, natomwfc
USE klist, ONLY: nks, xk
USE ldaU, ONLY: swfcatom
USE wvfct, ONLY: npwx, npw, igk
USE wvfct, ONLY: npwx, npw, igk, gamma_only
USE us, ONLY: nkb, vkb
use becmod
use becmod, only: becp
use rbecmod, only: rbecp => becp
#ifdef __PARA
use para
#endif
implicit none
!
IMPLICIT NONE
!
!
integer :: ik, ibnd, info, i, j, k, na, nb, nt, isym, n, ntemp, m, &
l, lm, ltot, ntot
integer, allocatable :: nwfc (:), ml (:)
@ -41,41 +45,52 @@ subroutine orthoatwfc
real(kind=DP) :: t0, scnds
! cpu time spent
logical :: orthogonalize_wfc
complex(kind=DP) :: ZDOTC, temp, t (5)
complex(kind=DP) , allocatable :: wfcatom (:,:), work (:,:), overlap (:,:)
real(kind=DP) , allocatable :: e (:)
t0 = scnds ()
allocate (wfcatom( npwx, natomwfc))
allocate (overlap( natomwfc , natomwfc))
allocate (work ( natomwfc , natomwfc))
allocate (e ( natomwfc))
allocate (ml ( natomwfc))
allocate (nwfc ( nat))
if ( gamma_only ) then
! Allocate the array rbecp=<beta|wfcatom>, which is not allocated by
! allocate_wfc() in a gamma-point calculation:
allocate (rbecp (nkb,natomwfc))
end if
orthogonalize_wfc = .false.
if (orthogonalize_wfc) then
WRITE( stdout,*) 'Atomic wfc used in LDA+U are orthogonalized'
WRITE( stdout,*) 'Gamma-only calculation for this case not implemented'
STOP
else
WRITE( stdout,*) 'Atomic wfc used in LDA+U are NOT orthogonalized'
end if
if (nks.gt.1) rewind (iunigk)
do ik = 1, nks
if (nks.gt.1) read (iunigk) npw, igk
overlap(:,:) = (0.d0,0.d0)
work(:,:) = (0.d0,0.d0)
call atomic_wfc (ik, wfcatom)
call init_us_2 (npw, igk, xk (1, ik), vkb)
call ccalbec (nkb, npwx, npw, natomwfc, becp, vkb, wfcatom)
if ( gamma_only ) then
call pw_gemm ('Y', nkb, natomwfc, npw, vkb, npwx, &
wfcatom, npwx, rbecp, nkb)
else
call ccalbec (nkb, npwx, npw, natomwfc, becp, vkb, wfcatom)
endif
call s_psi (npwx, npw, natomwfc, wfcatom, swfcatom)
@ -114,11 +129,11 @@ subroutine orthoatwfc
natomwfc, swfcatom (i, 1) , npwx, (0.d0, 0.d0) , work, 1)
call ZCOPY (natomwfc, work, 1, swfcatom (i, 1), npwx)
enddo
end if
call davcio (swfcatom, nwordatwfc, iunat, ik, 1)
enddo
deallocate (nwfc)
deallocate (ml)
@ -126,7 +141,12 @@ subroutine orthoatwfc
deallocate (work)
deallocate (e)
deallocate (wfcatom)
if ( gamma_only ) then
deallocate (rbecp)
end if
!
return
END SUBROUTINE orthoatwfc
end subroutine orthoatwfc

View File

@ -21,6 +21,8 @@ subroutine vhpsi (ldap, np, mp, psip, hpsi)
USE lsda_mod, ONLY: nspin, current_spin
USE basis, ONLY: nat, ntyp, ityp, natomwfc
USE us, ONLY: newpseudo
USE wvfct, ONLY: gamma_only
USE gvect, ONLY : gstart
implicit none
integer :: ldap, np, mp
@ -30,6 +32,7 @@ subroutine vhpsi (ldap, np, mp, psip, hpsi)
integer, allocatable :: offset (:)
! offset of localized electrons of atom na in the natomwfc ordering
complex(kind=DP) :: ZDOTC, temp
real(kind=DP), external :: DDOT
complex(kind=DP), allocatable :: proj (:,:)
!
allocate ( offset(nat), proj(natomwfc,mp) )
@ -48,7 +51,14 @@ subroutine vhpsi (ldap, np, mp, psip, hpsi)
if (counter.ne.natomwfc) call errore ('vhpsi', 'nstart<>counter', 1)
do ibnd = 1, mp
do i = 1, natomwfc
proj (i, ibnd) = ZDOTC (np, swfcatom (1, i), 1, psip (1, ibnd), 1)
if (gamma_only) then
proj (i, ibnd) = 2.d0 * &
DDOT(2*np, swfcatom (1, i), 1, psip (1, ibnd), 1)
if (gstart.eq.2) proj (i, ibnd) = proj (i, ibnd) - &
swfcatom (1, i) * psip (1, ibnd)
else
proj (i, ibnd) = ZDOTC (np, swfcatom (1, i), 1, psip (1, ibnd), 1)
endif
enddo
enddo
#ifdef __PARA
@ -67,12 +77,16 @@ subroutine vhpsi (ldap, np, mp, psip, hpsi)
temp = temp * Hubbard_U(nt)/2.d0
temp = temp + proj(offset(na)+m1,ibnd) * Hubbard_alpha(nt)
call ZAXPY (np, temp, swfcatom(1,offset(na)+m1), 1, &
if (gamma_only) then
call DAXPY (2*np, temp, swfcatom(1,offset(na)+m1), 1, &
hpsi(1,ibnd), 1)
else
call ZAXPY (np, temp, swfcatom(1,offset(na)+m1), 1, &
hpsi(1,ibnd), 1)
endif
enddo
endif
enddo
enddo
deallocate (offset, proj)
return

View File

@ -112,7 +112,7 @@ SUBROUTINE wfcinit()
!
! ... Needed for LDA+U (not yet in gamma)
!
!!!IF ( lda_plus_u ) CALL orthoatwfc()
IF ( lda_plus_u ) CALL orthoatwfc()
!
IF ( nks > 1 ) REWIND( iunigk )
!
@ -142,8 +142,8 @@ SUBROUTINE wfcinit()
!
g2kin(:) = g2kin(:) * tpiba2
!
!!!IF ( lda_plus_u ) &
!!! CALL davcio( swfcatom, nwordatwfc, iunat, ik, - 1 )
IF ( lda_plus_u ) &
CALL davcio( swfcatom, nwordatwfc, iunat, ik, - 1 )
!
IF ( startingwfc == 'atomic' ) THEN
!