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 .
|
|
|
|
!
|
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
|
|
|
#include "f_defs.h"
|
2004-12-21 23:48:19 +08:00
|
|
|
!-------------------------------------------------------------------------
|
|
|
|
subroutine calphiid(c0,bec,betae,phi)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
! input: c0 (orthonormal with s(r(t)), bec=<c0|beta>, betae=|beta>
|
|
|
|
! computes the matrix phi (with the old positions)
|
|
|
|
! where |phi> = s'|c0> = |c0> + sum q_ij |i><j|c0>
|
|
|
|
! where s'=s(r(t))
|
|
|
|
!
|
|
|
|
use ions_base, only: na, nsp
|
|
|
|
use io_global, only: stdout
|
|
|
|
use cvan
|
|
|
|
use uspp_param, only: nh
|
|
|
|
use uspp, only :nhsa=>nkb, nhsavb=>nkbus, qq
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: n => nbsp
|
2004-12-21 23:48:19 +08:00
|
|
|
use gvecw, only: ngw
|
|
|
|
use constants, only: pi, fpi
|
|
|
|
use control_flags, only: iprint, iprsta
|
|
|
|
!
|
|
|
|
implicit none
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) c0(ngw,n), phi(ngw,n), betae(ngw,nhsa)
|
|
|
|
real(8) bec(nhsa,n), emtot
|
2004-12-21 23:48:19 +08:00
|
|
|
! local variables
|
|
|
|
integer is, iv, jv, ia, inl, jnl, i, j
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) qtemp(nhsavb,n) ! automatic array
|
2004-12-21 23:48:19 +08:00
|
|
|
!
|
2005-10-22 03:40:15 +08:00
|
|
|
phi(1:ngw,1:n) = 0.0d0
|
2004-12-21 23:48:19 +08:00
|
|
|
!
|
|
|
|
if (nvb.gt.0) then
|
|
|
|
qtemp = 0.0d0
|
|
|
|
do is=1,nvb
|
|
|
|
do iv=1,nh(is)
|
|
|
|
do jv=1,nh(is)
|
|
|
|
if(abs(qq(iv,jv,is)).gt.1.e-5) then
|
|
|
|
do ia=1,na(is)
|
|
|
|
inl=ish(is)+(iv-1)*na(is)+ia
|
|
|
|
jnl=ish(is)+(jv-1)*na(is)+ia
|
|
|
|
do i=1,n
|
|
|
|
qtemp(inl,i) = qtemp(inl,i) + &
|
|
|
|
& qq(iv,jv,is)*bec(jnl,i)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
endif
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
!
|
|
|
|
call MXMA &
|
|
|
|
& (betae,1,2*ngw,qtemp,1,nhsavb,phi,1,2*ngw,2*ngw,nhsavb,n)
|
|
|
|
end if
|
|
|
|
!
|
|
|
|
do j=1,n
|
|
|
|
do i=1,ngw
|
|
|
|
phi(i,j)=(phi(i,j)+c0(i,j))
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
! =================================================================
|
|
|
|
if(iprsta.gt.2) then
|
|
|
|
emtot=0.
|
|
|
|
do j=1,n
|
|
|
|
do i=1,ngw
|
|
|
|
emtot=emtot &
|
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
|
|
|
& +2.*DBLE(phi(i,j)*CONJG(c0(i,j)))
|
2004-12-21 23:48:19 +08:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
emtot=emtot/n
|
|
|
|
#ifdef __PARA
|
|
|
|
call reduce(1,emtot)
|
|
|
|
#endif
|
|
|
|
WRITE( stdout,*) 'in calphi sqrt(emtot)=',sqrt(emtot)
|
|
|
|
WRITE( stdout,*)
|
|
|
|
do is=1,nsp
|
|
|
|
if(nsp.gt.1) then
|
|
|
|
WRITE( stdout,'(33x,a,i4)') ' calphi: bec (is)',is
|
|
|
|
WRITE( stdout,'(8f9.4)') &
|
|
|
|
& ((bec(ish(is)+(iv-1)*na(is)+1,i),iv=1,nh(is)),i=1,n)
|
|
|
|
else
|
|
|
|
do ia=1,na(is)
|
|
|
|
WRITE( stdout,'(33x,a,i4)') ' calphi: bec (ia)',ia
|
|
|
|
WRITE( stdout,'(8f9.4)') &
|
|
|
|
& ((bec(ish(is)+(iv-1)*na(is)+ia,i),iv=1,nh(is)),i=1,n)
|
|
|
|
end do
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
return
|
|
|
|
end subroutine calphiid
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine calcmt(fdiag,zmat,fmat)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! constructs fmat=zmat^t.fdiag.zmat
|
|
|
|
!
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: nudx, nspin, nupdwn, iupdwn, nx => nbspx
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer iss, nss, istart, i, j, k
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) zmat(nudx,nudx,nspin), fmat(nudx,nudx,nspin), &
|
2004-12-21 23:48:19 +08:00
|
|
|
& fdiag(nx)
|
|
|
|
|
|
|
|
do iss=1,nspin
|
|
|
|
nss=nupdwn(iss)
|
|
|
|
istart=iupdwn(iss)
|
|
|
|
do i=1,nss
|
|
|
|
do k=1,nss
|
|
|
|
fmat(k,i,iss)=0.0
|
|
|
|
do j=1,nss
|
|
|
|
fmat(k,i,iss)=fmat(k,i,iss)+ &
|
|
|
|
& zmat(j,k,iss)*fdiag(j+istart-1)*zmat(j,i,iss)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
return
|
2005-05-12 21:23:41 +08:00
|
|
|
end subroutine calcmt
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine rotate(z0,c0,bec,c0diag,becdiag)
|
|
|
|
!-----------------------------------------------------------------------
|
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
|
|
|
use kinds, only: dp
|
2004-12-21 23:48:19 +08:00
|
|
|
use cvan
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: nudx, nspin, nupdwn, iupdwn, nx => nbspx, n => nbsp
|
2004-12-21 23:48:19 +08:00
|
|
|
use uspp_param, only: nh
|
|
|
|
use uspp, only :nhsa=>nkb, nhsavb=>nkbus, qq
|
|
|
|
use gvecw, only: ngw
|
|
|
|
use ions_base, only: nas => nax, nsp, na
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer iss, nss, istart, i, j, k, ni, nj, is, jv, ia, jnl
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) z0(nudx,nudx,nspin)
|
|
|
|
real(8) bec(nhsa,n), becdiag(nhsa,n)
|
|
|
|
complex(8) c0(ngw,nx), c0diag(ngw,nx)
|
2004-12-21 23:48:19 +08:00
|
|
|
|
2005-10-22 03:40:15 +08:00
|
|
|
c0diag(1:ngw,1:nx)=0.d0
|
2004-12-21 23:48:19 +08:00
|
|
|
do iss=1,nspin
|
|
|
|
nss=nupdwn(iss)
|
|
|
|
istart=iupdwn(iss)
|
|
|
|
call MXMA(c0(1,istart),1,2*ngw,z0(1,1,iss),nudx,1,c0diag(1,istart),1,2*ngw,2*ngw,nss,nss)
|
|
|
|
end do
|
|
|
|
|
|
|
|
do iss=1,nspin
|
|
|
|
nss=nupdwn(iss)
|
|
|
|
istart=iupdwn(iss)
|
|
|
|
do ni=1,nss
|
|
|
|
do is=1,nsp
|
|
|
|
do jv=1,nh(is)
|
|
|
|
do ia=1,na(is)
|
|
|
|
jnl=ish(is)+(jv-1)*na(is)+ia
|
|
|
|
becdiag(jnl,ni+istart-1)=0.0
|
|
|
|
do nj=1,nss
|
|
|
|
becdiag(jnl,ni+istart-1)=becdiag(jnl,ni+istart-1)+ &
|
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
|
|
|
& CMPLX(z0(ni,nj,iss),0.d0)*bec(jnl,nj+istart-1)
|
2004-12-21 23:48:19 +08:00
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
return
|
2005-05-12 21:23:41 +08:00
|
|
|
end subroutine rotate
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine ddiag(nx,n,amat,dval,dvec,iflag)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
integer nx,n,naux,ndim,iopt,iflag,k,i,j,info
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) dval(n)
|
|
|
|
real(8) amat(nx,n), dvec(nx,n)
|
|
|
|
real(8), allocatable:: ap(:)
|
|
|
|
real(8), allocatable:: aux(:)
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
|
|
|
|
ndim=(n*(n+1))/2
|
|
|
|
naux=3*n
|
|
|
|
allocate(ap(ndim))
|
|
|
|
allocate(aux(naux))
|
|
|
|
|
|
|
|
|
|
|
|
if (iflag.eq.1) then
|
|
|
|
iopt=1
|
|
|
|
else if (iflag.eq.0) then
|
|
|
|
iopt=0
|
|
|
|
else
|
|
|
|
write(*,*) 'ERROR: diag, iflag not allowed'
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
|
|
|
|
k=0
|
|
|
|
do j=1,n
|
|
|
|
do i=1,j
|
|
|
|
k=k+1
|
|
|
|
! ap(i + (j-1)*j/2)=amat(i,j)
|
|
|
|
ap(k)=amat(i,j)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
call dspev('V','U',n,ap,dval,dvec,nx,aux,info)
|
2005-10-25 23:44:24 +08:00
|
|
|
if(info.ne.0) write(6,*) 'Problems with ddiag'
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
deallocate(ap)
|
|
|
|
deallocate(aux)
|
|
|
|
|
|
|
|
return
|
|
|
|
end subroutine ddiag
|
|
|
|
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine calcm(fdiag,zmat,fmat)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! constructs fmat=zmat.fdiag.zmat^t
|
|
|
|
!
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: nudx, nspin, nupdwn, iupdwn, nx => nbspx, n => nbsp
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer iss, nss, istart, i, j, k
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) zmat(nudx,nudx,nspin), fmat(nudx,nudx,nspin), &
|
2004-12-21 23:48:19 +08:00
|
|
|
& fdiag(nx)
|
|
|
|
|
|
|
|
do iss=1,nspin
|
|
|
|
nss=nupdwn(iss)
|
|
|
|
istart=iupdwn(iss)
|
|
|
|
do i=1,nss
|
|
|
|
do k=1,nss
|
|
|
|
fmat(k,i,iss)=0.0
|
|
|
|
do j=1,nss
|
|
|
|
fmat(k,i,iss)=fmat(k,i,iss)+ &
|
|
|
|
& zmat(k,j,iss)*fdiag(j+istart-1)*zmat(i,j,iss)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
return
|
2005-05-12 21:23:41 +08:00
|
|
|
end subroutine calcm
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
subroutine minparabola(ene0,dene0,ene1,passop,passo,stima)
|
2005-10-25 23:44:24 +08:00
|
|
|
!this subroutines finds the minimum of a quadratic real function
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
implicit none
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) ene0,dene0,ene1,passop,passo,stima
|
|
|
|
real(8) a,b,c!a*x^2+b*x+c
|
2004-12-21 23:48:19 +08:00
|
|
|
|
|
|
|
c=ene0
|
|
|
|
b=dene0
|
|
|
|
a=(ene1-b*passop-c)/(passop**2.d0)
|
|
|
|
|
|
|
|
passo = -b/(2.d0*a)
|
|
|
|
if( a.lt.0.d0) then
|
|
|
|
if(ene1.lt.ene0) then
|
|
|
|
passo=passop
|
|
|
|
else
|
|
|
|
passo=0.5*passop
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
|
|
|
stima=a*passo**2.d0+b*passo+c
|
|
|
|
|
|
|
|
|
|
|
|
return
|
|
|
|
end subroutine minparabola
|
|
|
|
|
|
|
|
|
|
|
|
subroutine pc2(a,beca,b,becb)
|
|
|
|
|
|
|
|
! this function applies the operator Pc
|
|
|
|
|
|
|
|
! this subroutine applies the Pc operator
|
|
|
|
! a input :unperturbed wavefunctions
|
|
|
|
! b input :first order wavefunctions
|
|
|
|
! b output:b_i =b_i-a_j><a_j|S|b_i>
|
|
|
|
|
|
|
|
use ions_base, only: na, nsp
|
|
|
|
use io_global, only: stdout
|
|
|
|
use cvan
|
|
|
|
use gvecw, only: ngw
|
|
|
|
use constants, only: pi, fpi
|
|
|
|
use control_flags, only: iprint, iprsta
|
|
|
|
use reciprocal_vectors, only: ng0 => gstart
|
|
|
|
use mp, only: mp_sum
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: n => nbsp
|
2004-12-21 23:48:19 +08:00
|
|
|
use uspp_param, only: nh
|
|
|
|
use uspp, only :nhsa=>nkb
|
|
|
|
use uspp, only :qq
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) a(ngw,n), b(ngw,n)
|
2004-12-21 23:48:19 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) beca(nhsa,n),becb(nhsa,n)
|
2004-12-21 23:48:19 +08:00
|
|
|
! local variables
|
|
|
|
integer is, iv, jv, ia, inl, jnl, i, j,ig
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) sca
|
2004-12-21 23:48:19 +08:00
|
|
|
do i=1,n
|
|
|
|
do j=1,n
|
|
|
|
sca=0.
|
|
|
|
if (ng0.eq.2) 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
|
|
|
b(1,i)=0.5d0*(b(1,i)+CONJG(b(1,i)))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
do ig=1,ngw !loop on g vectors
|
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
|
|
|
sca=sca+2.d0*DBLE(CONJG(a(ig,j))*b(ig,i)) !2. for real wavefunctions
|
2004-12-21 23:48:19 +08:00
|
|
|
enddo
|
|
|
|
if (ng0.eq.2) 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
|
|
|
sca=sca-DBLE(CONJG(a(1,j))*b(1,i))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
|
|
|
|
call mp_sum( sca )
|
|
|
|
|
|
|
|
if (nvb.gt.0) then
|
|
|
|
|
|
|
|
do is=1,nvb
|
|
|
|
do iv=1,nh(is)
|
|
|
|
do jv=1,nh(is)
|
|
|
|
do ia=1,na(is)
|
|
|
|
inl=ish(is)+(iv-1)*na(is)+ia
|
|
|
|
jnl=ish(is)+(jv-1)*na(is)+ia
|
|
|
|
sca=sca+ qq(iv,jv,is)*beca(inl,j)*becb(jnl,i)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
end if
|
|
|
|
|
|
|
|
do ig=1,ngw
|
|
|
|
b(ig,i)=b(ig,i)-sca*a(ig,j)
|
|
|
|
enddo
|
|
|
|
! this to prevent numerical errors
|
|
|
|
if (ng0.eq.2) 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
|
|
|
b(1,i)=0.5d0*(b(1,i)+CONJG(b(1,i)))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
|
|
|
return
|
|
|
|
end subroutine pc2
|
|
|
|
|
|
|
|
|
|
|
|
subroutine pcdaga2(a,as ,b )
|
|
|
|
|
|
|
|
! this function applies the operator Pc
|
|
|
|
|
|
|
|
! this subroutine applies the Pc^dagerr operator
|
|
|
|
! a input :unperturbed wavefunctions
|
|
|
|
! b input :first order wavefunctions
|
|
|
|
! b output:b_i =b_i - S|a_j><a_j|b_i>
|
|
|
|
|
|
|
|
use ions_base, only: na, nsp
|
|
|
|
use io_global, only: stdout
|
|
|
|
use cvan
|
|
|
|
use gvecw, only: ngw
|
|
|
|
use constants, only: pi, fpi
|
|
|
|
use control_flags, only: iprint, iprsta
|
|
|
|
use reciprocal_vectors, only: ng0 => gstart
|
|
|
|
use mp, only: mp_sum
|
2005-03-02 18:03:55 +08:00
|
|
|
use electrons_base, only: n => nbsp
|
2004-12-21 23:48:19 +08:00
|
|
|
use uspp_param, only: nh
|
|
|
|
use uspp, only :nhsa=>nkb
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(8) a(ngw,n), b(ngw,n), as(ngw,n)
|
|
|
|
real(8) beca(nhsa,n),becb(nhsa,n)
|
2004-12-21 23:48:19 +08:00
|
|
|
! local variables
|
|
|
|
integer is, iv, jv, ia, inl, jnl, i, j,ig
|
2005-08-28 22:09:42 +08:00
|
|
|
real(8) sca
|
2004-12-21 23:48:19 +08:00
|
|
|
!
|
|
|
|
do i=1,n
|
|
|
|
do j=1,n
|
|
|
|
sca=0.
|
|
|
|
if (ng0.eq.2) 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
|
|
|
b(1,i)=0.5d0*(b(1,i)+CONJG(b(1,i)))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
do ig=1,ngw !loop on g vectors
|
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
|
|
|
sca=sca+2.*DBLE(CONJG(a(ig,j))*b(ig,i)) !2. for real weavefunctions
|
2004-12-21 23:48:19 +08:00
|
|
|
enddo
|
|
|
|
if (ng0.eq.2) 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
|
|
|
sca=sca-DBLE(CONJG(a(1,j))*b(1,i))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
call mp_sum(sca)
|
|
|
|
do ig=1,ngw
|
|
|
|
b(ig,i)=b(ig,i)-sca*as(ig,j)
|
|
|
|
enddo
|
|
|
|
! this to prevent numerical errors
|
|
|
|
if (ng0.eq.2) 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
|
|
|
b(1,i)=0.5d0*(b(1,i)+CONJG(b(1,i)))
|
2004-12-21 23:48:19 +08:00
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
return
|
|
|
|
end subroutine pcdaga2
|
|
|
|
|