2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-09-27 17:11:56 +08:00
|
|
|
! Copyright (C) 2002-2004 PWSCF group
|
2003-01-20 05:58:50 +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-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
|
|
|
!-----------------------------------------------------------------------
|
2004-09-27 17:11:56 +08:00
|
|
|
SUBROUTINE dprojdtau(dproj,wfcatom,spsi,alpha,ipol,offset)
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! This routine computes the first derivative of the projection
|
2003-01-20 05:58:50 +08:00
|
|
|
! <\fi^{at}_{I,m1}|S|\psi_{k,v,s}> with respect to the atomic displacement
|
|
|
|
! u(alpha,ipol) (we remember that ns_{I,s,m1,m2} = \sum_{k,v}
|
|
|
|
! f_{kv} <\fi^{at}_{I,m1}|S|\psi_{k,v,s}><\psi_{k,v,s}|S|\fi^{at}_{I,m2}>)
|
|
|
|
!
|
2004-06-12 21:44:18 +08:00
|
|
|
USE kinds, ONLY : DP
|
|
|
|
USE atom, ONLY : nchi, lchi, oc
|
|
|
|
USE ions_base, ONLY : nat, ntyp => nsp, ityp
|
|
|
|
USE basis, ONLY : natomwfc
|
|
|
|
USE cell_base, ONLY : tpiba
|
|
|
|
USE gvect, ONLY : g
|
|
|
|
USE klist, ONLY : nks, xk
|
|
|
|
USE ldaU, ONLY : swfcatom, Hubbard_l, &
|
|
|
|
Hubbard_U, Hubbard_alpha
|
|
|
|
USE lsda_mod, ONLY : lsda, nspin, current_spin, isk
|
|
|
|
USE wvfct, ONLY : nbnd, npwx, npw, igk, wg
|
|
|
|
USE uspp, ONLY : nkb, vkb, qq
|
|
|
|
USE uspp_param, ONLY : nhm, nh
|
|
|
|
USE wavefunctions_module, ONLY : evc
|
2004-09-27 17:11:56 +08:00
|
|
|
USE becmod, ONLY : becp
|
|
|
|
|
|
|
|
IMPLICIT NONE
|
|
|
|
INTEGER :: &
|
2003-01-20 05:58:50 +08:00
|
|
|
alpha, &! input: the displaced atom
|
|
|
|
ipol, &! input: the component of displacement
|
|
|
|
offset ! input: the offset of the wfcs of the atom "alpha"
|
2004-09-27 17:11:56 +08:00
|
|
|
COMPLEX (kind=DP) :: &
|
2003-01-20 05:58:50 +08:00
|
|
|
wfcatom(npwx,natomwfc), &! input: the atomic wfc
|
|
|
|
spsi(npwx,nbnd), &! input: S|evc>
|
|
|
|
dproj(natomwfc,nbnd) ! output: the derivative of the projection
|
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
INTEGER :: ig, jkb2, na, m1, ibnd, iwf, nt, ib, ih,jh, ldim
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
REAL (kind=DP) :: gvec, a1, a2
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
COMPLEX (kind=DP):: ZDOTC
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
COMPLEX (kind=DP), ALLOCATABLE :: dwfc(:,:), work(:), dbeta(:), &
|
2003-01-20 05:58:50 +08:00
|
|
|
betapsi(:,:), dbetapsi(:,:), &
|
|
|
|
wfatbeta(:,:), wfatdbeta(:,:)
|
2003-02-10 16:58:33 +08:00
|
|
|
! dwfc(npwx,ldim), ! the derivative of the atomic d wfc
|
2003-01-20 05:58:50 +08:00
|
|
|
! work(npwx), ! the beta function
|
|
|
|
! dbeta(npwx), ! the derivative of the beta function
|
|
|
|
! betapsi(nhm,nbnd), ! <beta|evc>
|
|
|
|
! dbetapsi(nhm,nbnd), ! <dbeta|evc>
|
|
|
|
! wfatbeta(natomwfc,nhm),! <wfc|beta>
|
|
|
|
! wfatdbeta(natomwfc,nhm)! <wfc|dbeta>
|
|
|
|
|
2003-02-10 16:58:33 +08:00
|
|
|
nt = ityp(alpha)
|
|
|
|
|
|
|
|
ldim = 2 * Hubbard_l(nt) + 1
|
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
ALLOCATE ( dwfc(npwx,ldim), work(npwx), dbeta(npwx), betapsi(nhm,nbnd), &
|
2003-01-20 05:58:50 +08:00
|
|
|
dbetapsi(nhm,nbnd), wfatbeta(natomwfc,nhm), wfatdbeta(natomwfc,nhm) )
|
|
|
|
|
|
|
|
dproj(:,:) = (0.d0, 0.d0)
|
|
|
|
!
|
|
|
|
! At first the derivatives of the atomic wfc and the beta are computed
|
|
|
|
!
|
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
IF (Hubbard_U(nt).NE.0.d0.OR.Hubbard_alpha(nt).NE.0.d0) THEN
|
|
|
|
DO ig = 1,npw
|
2003-01-20 05:58:50 +08:00
|
|
|
gvec = g(ipol,igk(ig)) * tpiba
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2003-01-20 05:58:50 +08:00
|
|
|
! in the expression of dwfc we don't need (k+G) but just G; k always
|
|
|
|
! multiplies the underived quantity and gives an opposite contribution
|
|
|
|
! in c.c. term because the sign of the imaginary unit.
|
2003-02-10 16:58:33 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
DO m1 = 1, ldim
|
2003-01-20 05:58:50 +08:00
|
|
|
dwfc(ig,m1) = dcmplx(0.d0,-1.d0) * gvec * wfcatom(ig,offset+m1)
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL ZGEMM('C','N',ldim, nbnd, npw, (1.d0,0.d0), &
|
2003-02-10 16:58:33 +08:00
|
|
|
dwfc, npwx, spsi, npwx, (0.d0,0.d0), &
|
|
|
|
dproj(offset+1,1), natomwfc)
|
2004-09-27 17:11:56 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL reduce(2*natomwfc*nbnd,dproj)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
|
|
|
|
|
|
|
jkb2 = 0
|
2004-09-27 17:11:56 +08:00
|
|
|
DO nt=1,ntyp
|
|
|
|
DO na=1,nat
|
|
|
|
IF ( ityp(na) .EQ. nt ) THEN
|
|
|
|
DO ih=1,nh(nt)
|
2003-01-20 05:58:50 +08:00
|
|
|
jkb2 = jkb2 + 1
|
2004-09-27 17:11:56 +08:00
|
|
|
IF (na.EQ.alpha) THEN
|
|
|
|
DO ig = 1, npw
|
2003-01-20 05:58:50 +08:00
|
|
|
gvec = g(ipol,igk(ig)) * tpiba
|
2004-09-27 17:11:56 +08:00
|
|
|
dbeta(ig) = CMPLX(0.d0,-1.d0) * vkb(ig,jkb2) * gvec
|
2003-01-20 05:58:50 +08:00
|
|
|
work(ig) = vkb(ig,jkb2)
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
DO ibnd=1,nbnd
|
2003-01-20 05:58:50 +08:00
|
|
|
dbetapsi(ih,ibnd)= ZDOTC(npw,dbeta,1,evc(1,ibnd),1)
|
|
|
|
betapsi(ih,ibnd) = becp(jkb2,ibnd)
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
DO iwf=1,natomwfc
|
2003-01-20 05:58:50 +08:00
|
|
|
wfatbeta(iwf,ih) = ZDOTC(npw,wfcatom(1,iwf),1,work,1)
|
|
|
|
wfatdbeta(iwf,ih)= ZDOTC(npw,wfcatom(1,iwf),1,dbeta,1)
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
END IF
|
|
|
|
END DO
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL reduce(2*nhm*nbnd,dbetapsi)
|
|
|
|
CALL reduce(2*natomwfc*nhm,wfatbeta)
|
|
|
|
CALL reduce(2*natomwfc*nhm,wfatdbeta)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
2004-09-27 17:11:56 +08:00
|
|
|
IF (na.EQ.alpha) THEN
|
|
|
|
DO ibnd=1,nbnd
|
|
|
|
DO ih=1,nh(nt)
|
|
|
|
DO jh=1,nh(nt)
|
|
|
|
DO iwf=1,natomwfc
|
2003-01-20 05:58:50 +08:00
|
|
|
dproj(iwf,ibnd) = &
|
|
|
|
dproj(iwf,ibnd) + qq(ih,jh,nt) * &
|
|
|
|
( wfatdbeta(iwf,ih)*betapsi(jh,ibnd) + &
|
|
|
|
wfatbeta(iwf,ih)*dbetapsi(jh,ibnd) )
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END IF
|
|
|
|
END IF
|
|
|
|
END DO
|
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
DEALLOCATE ( dwfc, work, dbeta, betapsi, dbetapsi, wfatbeta, wfatdbeta )
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
RETURN
|
|
|
|
END SUBROUTINE dprojdtau
|