2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-01-28 20:28:11 +08:00
|
|
|
! Copyright (C) 2003 PWSCF group
|
|
|
|
! 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 .
|
|
|
|
!
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine A_h(e,h,ah)
|
|
|
|
!-----------------------------------------------------------------------
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only: DP
|
2004-05-12 05:08:21 +08:00
|
|
|
USE cell_base,ONLY : alat, omega, tpiba2
|
2004-06-01 01:55:33 +08:00
|
|
|
USE uspp, ONLY : vkb, nkb
|
2004-05-12 05:08:21 +08:00
|
|
|
USE lsda_mod, ONLY : current_spin, nspin
|
|
|
|
USE wvfct, ONLY: nbnd, npwx, npw, g2kin, igk
|
2003-11-09 18:42:50 +08:00
|
|
|
USE wavefunctions_module, ONLY: evc, psic
|
2004-05-12 05:08:21 +08:00
|
|
|
USE scf, ONLY : vrs, rho
|
|
|
|
USE gvect, ONLY : gstart, nl, nlm, nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
|
|
|
nrxx, ngm, g, gg
|
|
|
|
USE constants, ONLY: degspin, e2, fpi
|
2004-06-13 04:37:01 +08:00
|
|
|
use becmod, only: rbecp
|
2003-01-20 05:58:50 +08:00
|
|
|
use cgcom
|
|
|
|
use funct
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
integer :: j, jkb, ibnd, na,nt,ih
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: e(nbnd)
|
|
|
|
complex(DP) :: h(npwx,nbnd), ah(npwx,nbnd)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: fp, fm
|
|
|
|
complex(DP), pointer :: dpsic(:), drhoc(:), dvxc(:)
|
|
|
|
real(DP), pointer :: dv(:), drho(:)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
call start_clock('a_h')
|
|
|
|
!
|
|
|
|
drho => auxr
|
|
|
|
dpsic => aux2
|
|
|
|
drhoc => aux3
|
|
|
|
!
|
2003-09-18 05:50:03 +08:00
|
|
|
drho(:) = 0.d0
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! [(k+G)^2 - e ]psi
|
|
|
|
do ibnd = 1,nbnd
|
2003-02-08 00:04:36 +08:00
|
|
|
! set to zero the imaginary part of h at G=0
|
2003-01-20 05:58:50 +08:00
|
|
|
! needed for numerical stability
|
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
|
|
|
if (gstart==2) h(1,ibnd) = CMPLX( DBLE(h(1,ibnd)),0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
do j = 1,npw
|
|
|
|
ah(j,ibnd) = (g2kin(j)-e(ibnd)) * h(j,ibnd)
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
! V_Loc psi
|
|
|
|
do ibnd = 1,nbnd, 2
|
2003-09-18 05:50:03 +08:00
|
|
|
dpsic(:)= (0.d0, 0.d0)
|
|
|
|
psic(:) = (0.d0, 0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
if (ibnd.lt.nbnd) then
|
|
|
|
! two ffts at the same time
|
|
|
|
do j = 1,npw
|
|
|
|
psic (nl (igk(j))) = evc(j,ibnd) + (0.0,1.d0)* evc(j,ibnd+1)
|
|
|
|
dpsic(nl (igk(j))) = h(j,ibnd) + (0.0,1.d0)* h(j,ibnd+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
|
|
|
psic (nlm(igk(j)))= CONJG(evc(j,ibnd)-(0.0,1.d0)* evc(j,ibnd+1))
|
|
|
|
dpsic(nlm(igk(j)))= CONJG( h(j,ibnd)-(0.0,1.d0)* h(j,ibnd+1))
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
else
|
|
|
|
do j = 1,npw
|
|
|
|
psic (nl (igk(j))) = evc(j,ibnd)
|
|
|
|
dpsic(nl (igk(j))) = h(j,ibnd)
|
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
|
|
|
psic (nlm(igk(j))) = CONJG( evc(j,ibnd))
|
|
|
|
dpsic(nlm(igk(j))) = CONJG( h(j,ibnd))
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
end if
|
|
|
|
call cft3s( psic,nr1,nr2,nr3,nrx1,nr2,nr3,2)
|
|
|
|
call cft3s(dpsic,nr1,nr2,nr3,nrx1,nr2,nr3,2)
|
|
|
|
do j = 1,nrxx
|
|
|
|
drho(j) = drho(j) - 2.0*degspin/omega * &
|
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(psic(j)*CONJG(dpsic(j)))
|
2003-01-20 05:58:50 +08:00
|
|
|
dpsic(j) = dpsic(j) * vrs(j,current_spin)
|
|
|
|
end do
|
|
|
|
call cft3s(dpsic,nr1,nr2,nr3,nrx1,nr2,nr3,-2)
|
|
|
|
if (ibnd.lt.nbnd) then
|
|
|
|
! two ffts at the same time
|
|
|
|
do j = 1,npw
|
|
|
|
fp = (dpsic (nl(igk(j))) + dpsic (nlm(igk(j))))*0.5d0
|
|
|
|
fm = (dpsic (nl(igk(j))) - dpsic (nlm(igk(j))))*0.5d0
|
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
|
|
|
ah(j,ibnd ) = ah(j,ibnd) +CMPLX( DBLE(fp), AIMAG(fm))
|
|
|
|
ah(j,ibnd+1) = ah(j,ibnd+1)+CMPLX(AIMAG(fp),- DBLE(fm))
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
else
|
|
|
|
do j = 1,npw
|
|
|
|
ah(j,ibnd) = ah(j,ibnd) + dpsic (nl(igk(j)))
|
|
|
|
end do
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
!
|
|
|
|
nullify(dpsic)
|
|
|
|
! V_NL psi
|
2004-06-13 04:37:01 +08:00
|
|
|
call pw_gemm ('Y', nkb, nbnd, npw, vkb, npwx, h, npwx, rbecp, nkb)
|
2003-01-20 05:58:50 +08:00
|
|
|
if (nkb.gt.0) call add_vuspsi (npwx, npw, nbnd, h, ah)
|
|
|
|
!
|
|
|
|
do j = 1,nrxx
|
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
|
|
|
drhoc(j) = CMPLX(drho(j),0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
call cft3(drhoc,nr1,nr2,nr3,nrx1,nr2,nr3,-1)
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! drho is deltarho(r), drhoc is deltarho(g)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! mu'(n(r)) psi(r) delta psi(r)
|
|
|
|
!
|
|
|
|
dvxc => aux2
|
|
|
|
do j = 1,nrxx
|
|
|
|
dvxc(j) = drho(j)*dmuxc(j)
|
|
|
|
end do
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! add gradient correction contribution (if any)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
call start_clock('dgradcorr')
|
|
|
|
if (igcx.ne.0.or.igcc.ne.0) call dgradcor1 &
|
|
|
|
(rho, grho, dvxc_rr, dvxc_sr, dvxc_ss, dvxc_s, &
|
|
|
|
drho, drhoc, nr1,nr2,nr3, nrx1, nrx2, nrx3, nrxx, nspin, &
|
|
|
|
nl, nlm, ngm, g, alat, omega, dvxc)
|
|
|
|
call stop_clock('dgradcorr')
|
|
|
|
nullify (drho)
|
|
|
|
!
|
|
|
|
! 1/|r-r'| * psi(r') delta psi(r')
|
|
|
|
!
|
|
|
|
! gstart is the first nonzero G vector (needed for parallel execution)
|
|
|
|
!
|
|
|
|
if (gstart==2) drhoc(nl(1)) = 0.d0
|
|
|
|
!
|
|
|
|
do j = gstart,ngm
|
|
|
|
drhoc(nl (j)) = e2*fpi*drhoc(nl(j))/ (tpiba2*gg(j))
|
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
|
|
|
drhoc(nlm(j)) = CONJG(drhoc(nl (j)))
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
call cft3(drhoc,nr1,nr2,nr3,nrx1,nr2,nr3,+1)
|
|
|
|
!
|
|
|
|
! drhoc now contains deltaV_hartree
|
|
|
|
!
|
|
|
|
dv => auxr
|
|
|
|
do j = 1,nrxx
|
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
|
|
|
dv(j) = - DBLE(dvxc(j)) - DBLE(drhoc(j))
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call vloc_psi(npwx, npw, nbnd, evc, dv, ah)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! set to zero the imaginary part of ah at G=0
|
|
|
|
! needed for numerical stability
|
|
|
|
if (gstart.eq.2) then
|
|
|
|
do ibnd = 1, nbnd
|
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
|
|
|
ah(1,ibnd) = CMPLX( DBLE(ah(1,ibnd)),0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
end do
|
|
|
|
end if
|
|
|
|
!
|
|
|
|
call stop_clock('a_h')
|
|
|
|
!
|
|
|
|
return
|
|
|
|
end subroutine A_h
|