quantum-espresso/Gamma/h_h.f90

58 lines
1.5 KiB
Fortran
Raw Normal View History

!
! 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 .
!
!
!-----------------------------------------------------------------------
subroutine H_h(e,h,Ah)
!-----------------------------------------------------------------------
#include "f_defs.h"
!
USE kinds, only: DP
USE wvfct, ONLY: nbnd, npwx, npw, g2kin, igk
USE gvect, ONLY : gstart
USE uspp, ONLY : vkb, nkb
USE lsda_mod, ONLY : current_spin
USE scf, ONLY : vrs
use becmod, only: rbecp
use cgcom
!
implicit none
!
real(DP):: e(nbnd)
complex(DP):: h(npwx,nbnd), Ah(npwx,nbnd)
!
integer:: j,ibnd
!
call start_clock('h_h')
!
! [(k+G)^2 - e ]psi
do ibnd = 1,nbnd
! set to zero the imaginary part of h at G=0
! 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)
do j = 1,npw
ah(j,ibnd) = (g2kin(j)-e(ibnd)) * h(j,ibnd)
end do
end do
! V_Loc psi
call vloc_psi(npwx, npw, nbnd, h, vrs(1,current_spin), ah)
! V_NL psi
call pw_gemm ('Y', nkb, nbnd, npw, vkb, npwx, h, npwx, rbecp, nkb)
if (nkb.gt.0) call add_vuspsi (npwx, npw, nbnd, h, ah)
! set to zero the imaginary part of ah at G=0
! needed for numerical stability
if (gstart==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)
end do
end if
!
call stop_clock('h_h')
!
return
end subroutine H_h