quantum-espresso/PH/compute_becalp.f90

72 lines
2.0 KiB
Fortran
Raw Normal View History

!
! Copyright (C) 2001 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 compute_becalp (becq, alpq)
!---------------------------------------------------------------------
!
! This routine is used only if .not.lgamma and in this case
! computes the scalar product of vkb and psi_{k+q}
!
#include "f_defs.h"
use pwcom
USE kinds, only : DP
USE io_files, ONLY: iunigk
use phcom
implicit none
complex(DP) :: becq(nkb, nbnd, nksq), alpq(nkb,nbnd,3,nksq)
! the becp with psi_{k+q}
! the alphap with psi_{k+q}
integer :: ik, ikq, ipol, ibnd, ig, ios
! counter on k points
! counter on polarizations, bands and
! used for i/o control
complex(DP) :: fact
complex(DP), allocatable :: aux (:,:)
!
if (lgamma) return
allocate (aux ( npwx , nbnd))
if (nksq.gt.1) rewind (iunigk)
do ik = 1, nksq
ikq = 2 * ik
if (nksq.gt.1) then
read (iunigk, err = 100, iostat = ios) npw, igk
100 call errore ('compute_becalp', 'reading igk', abs (ios) )
read (iunigk, err = 200, iostat = ios) npwq, igkq
200 call errore ('compute_becalp', 'reading igkq', abs (ios) )
endif
call init_us_2 (npwq, igkq, xk (1, ikq), vkb)
call davcio (evq, lrwfc, iuwfc, ikq, - 1)
call ccalbec (nkb, npwx, npwq, nbnd, becq (1, 1, ik), vkb, evq)
do ipol = 1, 3
do ibnd = 1, nbnd
do ig = 1, npwq
aux (ig, ibnd) = evq (ig, ibnd) * &
(xk (ipol, ikq) + g (ipol, igkq(ig) ) )
enddo
enddo
call ccalbec (nkb, npwx, npwq, nbnd, alpq(1, 1, ipol, ik),vkb, aux)
enddo
enddo
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
fact = CMPLX (0.d0, tpiba)
call ZSCAL (nkb * nbnd * 3 * nksq, fact, alpq, 1)
deallocate (aux)
return
end subroutine compute_becalp