2005-05-17 03:19:04 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2002-2005 FPMD-CPV groups
|
|
|
|
! 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 .
|
|
|
|
!
|
|
|
|
subroutine dforceb(c0, i, betae, ipol, bec0, ctabin, gqq, gqqm, qmat, dq2, df)
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
! this subroutine computes the force for electrons
|
|
|
|
! in case of Berry,s phase like perturbation
|
|
|
|
! see internal notes
|
|
|
|
! it gives the force for the i-th state
|
|
|
|
!respect to vectorial (serial) program I changed ngwx to ngw :-)
|
|
|
|
|
|
|
|
! c0 input: Psi^0_i
|
|
|
|
! c1 input: Psi^1_i
|
|
|
|
! i input: ot computes the force for the i-th state
|
|
|
|
! v0 input: the local zeroth order potential
|
|
|
|
! v1 input: the local first order potential
|
|
|
|
! betae input: the functions beta_iR
|
|
|
|
! ipol input:the polarization of nabla_k
|
|
|
|
! bec0 input: the factors <beta_iR|Psi^0_v>
|
|
|
|
! bec1 input: the factors <beta_iR|Psi^1_v>
|
|
|
|
! ctabin input: the inverse-correspondence array g'+(-)1=g
|
|
|
|
! gqq input: the factors int dr Beta_Rj*Beta_Ri exp(iGr)
|
|
|
|
! gqqm input: the factors int dr Beta_Rj*Beta_Ri exp(iGr)
|
|
|
|
! qmat input:
|
|
|
|
! dq2 input: factors d^2hxc_ijR
|
|
|
|
! df output: force for the i-th state
|
|
|
|
|
|
|
|
|
|
|
|
use gvecs
|
|
|
|
use gvecw, only: ngw
|
|
|
|
use parameters
|
2005-04-15 05:08:53 +08:00
|
|
|
use electrons_base, only: nx => nbspx, n => nbsp, nspin, f
|
2004-12-21 23:48:19 +08:00
|
|
|
use constants
|
|
|
|
use cvan
|
|
|
|
use ions_base
|
|
|
|
use ions_base, only : nas => nax
|
|
|
|
use cell_base, only: a1, a2, a3
|
|
|
|
use uspp_param, only: nh, nhm
|
|
|
|
use uspp, only : nhsa=> nkb
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) c0(ngw, n), betae(ngw,nhsa), df(ngw),&
|
2004-12-21 23:48:19 +08:00
|
|
|
& gqq(nhm,nhm,nas,nsp),gqqm(nhm,nhm,nas,nsp),&
|
|
|
|
& qmat(nx,nx)
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) bec0(nhsa,n),&
|
2004-12-21 23:48:19 +08:00
|
|
|
& dq2(nat,nhm,nhm,nspin), gmes
|
|
|
|
|
|
|
|
integer i, ipol, ctabin(ngw,2)
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
|
|
|
|
integer j,k,ig,iv,jv,ix,jx,is,ia, isa,iss,iss1,mism
|
|
|
|
integer ir,ism,itemp,itempa,jnl,inl
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) ci ,fi, fp, fm
|
|
|
|
real(8) afr(nhsa), dd
|
|
|
|
complex(8) afrc(nhsa)
|
|
|
|
complex(8), allocatable:: dtemp(:)
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
allocate( dtemp(ngw))
|
|
|
|
|
|
|
|
|
|
|
|
ci=(0.,1.)
|
|
|
|
|
|
|
|
|
|
|
|
! now the interaction term
|
|
|
|
! first the norm-conserving part
|
|
|
|
|
|
|
|
do ig=1,ngw
|
|
|
|
dtemp(ig)=(0.,0.)
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do j=1,n
|
|
|
|
do ig=1,ngw
|
|
|
|
if(ctabin(ig,2) .ne. (ngw+1)) then
|
|
|
|
if(ctabin(ig,2).ge.0) then
|
|
|
|
dtemp(ig)=dtemp(ig)+c0(ctabin(ig,2),j)*qmat(j,i)
|
|
|
|
else
|
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
|
|
|
dtemp(ig)=dtemp(ig)+CONJG(c0(-ctabin(ig,2),j))*qmat(j,i)
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
do ig=1,ngw
|
|
|
|
if(ctabin(ig,1) .ne. (ngw+1)) then
|
|
|
|
if(ctabin(ig,1).ge.0) then
|
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
|
|
|
dtemp(ig)=dtemp(ig)-c0(ctabin(ig,1),j)*CONJG(qmat(j,i))
|
2004-12-21 23:48:19 +08:00
|
|
|
else
|
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
|
|
|
dtemp(ig)=dtemp(ig)-CONJG(c0(-ctabin(ig,1),j))*conjg(qmat(j,i))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
|
|
|
|
if(ipol.eq.1) then
|
|
|
|
gmes=a1(1)**2+a1(2)**2+a1(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
if(ipol.eq.2) then
|
|
|
|
gmes=a2(1)**2+a2(2)**2+a2(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
if(ipol.eq.3) then
|
|
|
|
gmes=a3(1)**2+a3(2)**2+a3(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
|
|
|
|
fi=f(i)*ci/(2.*gmes)
|
|
|
|
|
|
|
|
do ig=1,ngw
|
|
|
|
df(ig)= fi*dtemp(ig)
|
|
|
|
end do
|
|
|
|
|
|
|
|
! now the interacting Vanderbilt term
|
|
|
|
! the term (-ie/|G|)(-beta_i'R>gqq(i',j')bec0_jRj'Q^-1_ji+
|
|
|
|
! +beta_i'R>gqqm(i',j')bec0jRj'Q^-1_ij*
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
if(nhsa.gt.0) then
|
|
|
|
do inl=1,nhsa
|
|
|
|
afrc(inl)=(0.,0.)
|
|
|
|
end do
|
|
|
|
|
|
|
|
do is=1,nvb!loop on species
|
|
|
|
do iv=1,nh(is) !loop on projectors
|
|
|
|
do jv=1,nh(is) !loop on projectors
|
|
|
|
do ia=1,na(is)
|
|
|
|
inl=ish(is)+(iv-1)*na(is)+ia
|
|
|
|
jnl=ish(is)+(jv-1)*na(is)+ia
|
|
|
|
do j=1,n !loop on states
|
|
|
|
afrc(inl)=afrc(inl)+gqq(iv,jv,ia,is)*bec0(jnl,j)*qmat(j,i)&
|
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
|
|
|
& -CONJG(gqq(jv,iv,ia,is))*bec0(jnl,j)*conjg(qmat(i,j))
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do ig=1,ngw
|
|
|
|
dtemp(ig)=(0.,0.)
|
|
|
|
end do
|
|
|
|
do inl=1,nhsa
|
|
|
|
do ig=1,ngw
|
|
|
|
dtemp(ig)=dtemp(ig)+afrc(inl)*betae(ig,inl)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
! call MXMA
|
|
|
|
! & (betae,1,2*ngw,afr,1,nhsax,dtemp,1,2*ngw,2*ngw,nhsa,1)
|
|
|
|
do ig=1,ngw
|
|
|
|
df(ig)=df(ig)+fi*dtemp(ig)
|
|
|
|
end do
|
|
|
|
endif
|
|
|
|
|
|
|
|
deallocate( dtemp)
|
|
|
|
return
|
|
|
|
end subroutine dforceb
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
subroutine enberry( detq, ipol, enb)
|
|
|
|
|
|
|
|
use constants
|
|
|
|
use parameters
|
|
|
|
use cell_base, only: a1, a2, a3
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) detq
|
|
|
|
real(8) enb
|
2004-12-21 23:48:19 +08:00
|
|
|
integer ipol
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) gmes
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
|
|
|
|
if(ipol.eq.1) then
|
|
|
|
gmes=a1(1)**2+a1(2)**2+a1(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
if(ipol.eq.2) then
|
|
|
|
gmes=a2(1)**2+a2(2)**2+a2(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
if(ipol.eq.3) then
|
|
|
|
gmes=a3(1)**2+a3(2)**2+a3(3)**2
|
|
|
|
gmes=2*pi/SQRT(gmes)
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
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
|
|
|
enb = 2.*AIMAG(log(detq))/gmes!ATTENZIONE al segno
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
! write(6,*) detq, enb
|
|
|
|
return
|
|
|
|
end subroutine enberry
|
|
|
|
|
|
|
|
|
|
|
|
|