2003-01-20 05:58:50 +08:00
|
|
|
!
|
2007-10-24 03:47:26 +08:00
|
|
|
! Copyright (C) 2002-2007 Quantum-Espresso 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 dndepsilon ( dns,ldim,ipol,jpol )
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
! This routine computes the derivative of the ns atomic occupations with
|
|
|
|
! respect to the strain epsilon(ipol,jpol) used to obtain the hubbard
|
|
|
|
! contribution to the internal stres tensor.
|
|
|
|
!
|
2004-06-12 21:44:18 +08:00
|
|
|
USE kinds, ONLY : DP
|
|
|
|
USE wavefunctions_module, ONLY : evc
|
|
|
|
USE ions_base, ONLY : nat, ityp
|
|
|
|
USE basis, ONLY : natomwfc
|
2007-01-23 00:38:47 +08:00
|
|
|
USE klist, ONLY : nks, xk, ngk
|
2004-06-12 21:44:18 +08:00
|
|
|
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
|
2007-10-24 03:47:26 +08:00
|
|
|
USE uspp_param, ONLY : upf
|
2007-12-05 22:14:12 +08:00
|
|
|
USE becmod, ONLY : becp, calbec
|
2006-11-05 11:02:28 +08:00
|
|
|
USE io_files, ONLY : iunigk, nwordwfc, iunwfc, &
|
|
|
|
iunat, iunsat, nwordatwfc
|
2007-02-21 21:01:31 +08:00
|
|
|
USE buffers, ONLY : get_buffer
|
2008-01-24 00:53:17 +08:00
|
|
|
USE mp_global, ONLY : intra_pool_comm, inter_pool_comm
|
|
|
|
USE mp, ONLY : mp_sum
|
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
IMPLICIT NONE
|
2003-02-08 00:04:36 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
! I/O variables first
|
|
|
|
!
|
2004-09-27 17:11:56 +08:00
|
|
|
INTEGER :: ipol, jpol, ldim
|
2005-08-28 22:09:42 +08:00
|
|
|
REAL (DP) :: dns(ldim,ldim,nspin,nat)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! local variable
|
|
|
|
!
|
2004-09-27 17:11:56 +08:00
|
|
|
INTEGER :: ik, & ! counter on k points
|
2003-01-20 05:58:50 +08:00
|
|
|
ibnd, & ! " " bands
|
|
|
|
is, & ! " " spins
|
|
|
|
i, na, nt, n, counter, m1, m2, l
|
2005-08-28 22:09:42 +08:00
|
|
|
COMPLEX (DP) :: ZDOTC
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
INTEGER, ALLOCATABLE :: offset(:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! offset(nat) ! offset of d electrons of atom d in the natomwfc ordering
|
2005-08-28 22:09:42 +08:00
|
|
|
COMPLEX (DP), ALLOCATABLE :: &
|
2003-01-20 05:58:50 +08:00
|
|
|
proj(:,:), wfcatom(:,:), spsi(:,:), dproj(:,:)
|
2003-02-08 00:04:36 +08:00
|
|
|
! proj(natomwfc,nbnd), wfcatom(npwx,natomwfc),
|
2003-01-20 05:58:50 +08:00
|
|
|
! spsi(npwx,nbnd), dproj(natomwfc,nbnd)
|
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
ALLOCATE (offset(nat), proj(natomwfc,nbnd), wfcatom(npwx,natomwfc), &
|
2004-05-12 05:08:21 +08:00
|
|
|
spsi(npwx,nbnd), dproj(natomwfc,nbnd), becp(nkb,nbnd) )
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! D_Sl for l=1 and l=2 are already initialized, for l=0 D_S0 is 1
|
|
|
|
!
|
|
|
|
counter = 0
|
2004-09-27 17:11:56 +08:00
|
|
|
DO na=1,nat
|
2003-01-20 05:58:50 +08:00
|
|
|
offset(na) = 0
|
|
|
|
nt=ityp(na)
|
2007-10-24 03:47:26 +08:00
|
|
|
DO n=1,upf(nt)%nwfc
|
|
|
|
IF (upf(nt)%oc(n) >= 0.d0) THEN
|
|
|
|
l=upf(nt)%lchi(n)
|
2004-09-27 17:11:56 +08:00
|
|
|
IF (l.EQ.Hubbard_l(nt)) offset(na) = counter
|
2003-01-20 05:58:50 +08:00
|
|
|
counter = counter + 2 * l + 1
|
2004-09-27 17:11:56 +08:00
|
|
|
END IF
|
|
|
|
END DO
|
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
IF(counter.NE.natomwfc) CALL errore('new_ns','nstart<>counter',1)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2003-01-20 05:58:50 +08:00
|
|
|
dns(:,:,:,:) = 0.d0
|
|
|
|
!
|
|
|
|
! we start a loop on k points
|
|
|
|
!
|
2007-01-23 00:38:47 +08:00
|
|
|
IF (nks > 1) REWIND (iunigk)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
DO ik = 1, nks
|
|
|
|
IF (lsda) current_spin = isk(ik)
|
2007-01-23 00:38:47 +08:00
|
|
|
IF (nks > 1) READ (iunigk) igk
|
|
|
|
npw = ngk(ik)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! now we need the first derivative of proj with respect to
|
2003-01-20 05:58:50 +08:00
|
|
|
! epsilon(ipol,jpol)
|
|
|
|
!
|
2007-02-21 21:01:31 +08:00
|
|
|
CALL get_buffer (evc, nwordwfc, iunwfc, ik)
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL init_us_2 (npw,igk,xk(1,ik),vkb)
|
2007-12-05 22:14:12 +08:00
|
|
|
CALL calbec( npw, vkb, evc, becp )
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL s_psi (npwx, npw, nbnd, evc, spsi )
|
2006-11-05 11:02:28 +08:00
|
|
|
! CALL atomic_wfc( ik, wfcatom )
|
|
|
|
! read atomic wfc instead
|
|
|
|
CALL davcio(wfcatom,nwordatwfc,iunat,ik,-1)
|
|
|
|
|
2003-01-20 05:58:50 +08:00
|
|
|
dproj(:,:) = (0.d0,0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
CALL dprojdepsilon(ik,dproj,wfcatom,spsi,ipol,jpol)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2006-11-05 10:47:07 +08:00
|
|
|
CALL davcio(swfcatom,nwordatwfc,iunsat,ik,-1)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
DO ibnd = 1, nbnd
|
|
|
|
DO i=1,natomwfc
|
2003-01-20 05:58:50 +08:00
|
|
|
proj(i,ibnd) = ZDOTC(npw,swfcatom(1,i),1,evc(1,ibnd),1)
|
2004-09-27 17:11:56 +08:00
|
|
|
ENDDO
|
|
|
|
ENDDO
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2008-04-21 05:23:37 +08:00
|
|
|
CALL mp_sum( proj, intra_pool_comm )
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! compute the derivative of the occupation numbers (quantities dn(m1,m2))
|
|
|
|
! of the atomic orbitals. They are real quantities as well as n(m1,m2)
|
|
|
|
!
|
2004-09-27 17:11:56 +08:00
|
|
|
DO na = 1,nat
|
2003-01-20 05:58:50 +08:00
|
|
|
nt = ityp(na)
|
2004-09-27 17:11:56 +08:00
|
|
|
IF (Hubbard_U(nt).NE.0.d0.OR.Hubbard_alpha(nt).NE.0.d0) THEN
|
|
|
|
DO m1 = 1, 2 * Hubbard_l(nt) + 1
|
|
|
|
DO m2 = m1, 2 * Hubbard_l(nt) + 1
|
|
|
|
DO ibnd = 1,nbnd
|
2003-07-31 20:54:09 +08:00
|
|
|
dns(m1,m2,current_spin,na) = dns(m1,m2,current_spin,na) + &
|
2003-01-20 05:58:50 +08:00
|
|
|
wg(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
|
|
|
DBLE( proj(offset(na)+m1,ibnd) * &
|
2004-09-27 17:11:56 +08:00
|
|
|
CONJG(dproj(offset(na)+m2,ibnd) ) + &
|
2003-01-20 05:58:50 +08:00
|
|
|
dproj(offset(na)+m1,ibnd)* &
|
2004-09-27 17:11:56 +08:00
|
|
|
CONJG( proj(offset(na)+m2,ibnd) ) )
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END IF
|
|
|
|
END DO
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO ! on k-points
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2008-01-24 00:53:17 +08:00
|
|
|
CALL mp_sum( dns, inter_pool_comm )
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
2004-12-02 23:53:17 +08:00
|
|
|
!
|
|
|
|
! In nspin.eq.1 k-point weight wg is normalized to 2 el/band
|
|
|
|
! in the whole BZ but we are interested in dns of one spin component
|
|
|
|
!
|
|
|
|
IF (nspin.EQ.1) dns = 0.5d0 * dns
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! impose hermeticity of dn_{m1,m2}
|
|
|
|
!
|
2004-09-27 17:11:56 +08:00
|
|
|
DO na = 1,nat
|
2003-02-10 16:58:33 +08:00
|
|
|
nt = ityp(na)
|
2004-09-27 17:11:56 +08:00
|
|
|
DO is = 1,nspin
|
|
|
|
DO m1 = 1, 2 * Hubbard_l(nt) + 1
|
|
|
|
DO m2 = m1+1, 2 * Hubbard_l(nt) + 1
|
2003-07-31 20:54:09 +08:00
|
|
|
dns(m2,m1,is,na) = dns(m1,m2,is,na)
|
2004-09-27 17:11:56 +08:00
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END DO
|
|
|
|
END DO
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
DEALLOCATE (offset, proj, wfcatom, spsi, dproj, becp )
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-09-27 17:11:56 +08:00
|
|
|
RETURN
|
|
|
|
END SUBROUTINE dndepsilon
|