quantum-espresso/PH/addusdbec.f90

120 lines
3.3 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 addusdbec (ik, wgt, psi, dbecsum)
!----------------------------------------------------------------------
!
! This routine adds to the dbecsum the term which correspond to this
! k point. After the accumulation the additional part of the charge
! is computed in addusddens.
!
#include "f_defs.h"
USE ions_base, ONLY : nat, ityp, ntyp => nsp
use pwcom
USE kinds, only : DP
USE uspp_param, only: nh, tvanp, nhm
use phcom
implicit none
!
! the dummy variables
!
complex(DP) :: dbecsum (nhm*(nhm+1)/2, nat), psi(npwx,nbnd)
! inp/out: the sum kv of bec *
! input : contains delta psi
integer :: ik
! input: the k point
real(DP) :: wgt
! input: the weight of this k point
!
! here the local variables
!
integer :: na, nt, ih, jh, ibnd, ikk, ikb, jkb, ijh, startb, &
lastb, ijkb0
! counter on atoms
! counter on atomic type
! counter on solid beta functions
! counter on solid beta functions
! counter on the bands
! the real k point
! counter on solid becp
! counter on solid becp
! composite index for dbecsum
! divide among processors the sum
! auxiliary variable for counting
real(DP) :: w
! the weight of the k point
complex(DP), allocatable :: dbecq (:,:)
! the change of becq
if (.not.okvan) return
call start_clock ('addusdbec')
allocate (dbecq( nkb, nbndx))
w = 0.5d0 * wgt * omega
if (lgamma) then
ikk = ik
else
ikk = 2 * ik - 1
endif
!
! First compute the product of psi and vkb
!
call ccalbec (nkb, npwx, npwq, nbnd, dbecq, vkb, psi)
!
! And then we add the product to becsum
!
! Band parallelization: each processor takes care of its slice of bands
!
call divide (nbnd_occ (ikk), startb, lastb)
!
ijkb0 = 0
do nt = 1, ntyp
if (tvanp (nt) ) then
do na = 1, nat
if (ityp (na) .eq.nt) then
!
! And qgmq and becp and dbecq
!
ijh = 1
do ih = 1, nh (nt)
ikb = ijkb0 + ih
do ibnd = startb, lastb
dbecsum (ijh, na) = dbecsum (ijh, na) + &
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
w * ( CONJG(becp1(ikb,ibnd,ik)) * dbecq(ikb,ibnd) )
enddo
ijh = ijh + 1
do jh = ih + 1, nh (nt)
jkb = ijkb0 + jh
do ibnd = startb, lastb
dbecsum (ijh, na) = dbecsum (ijh, na) + &
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
w * ( CONJG(becp1(ikb,ibnd,ik)) * dbecq(jkb,ibnd) + &
CONJG(becp1(jkb,ibnd,ik)) * dbecq(ikb,ibnd) )
enddo
ijh = ijh + 1
enddo
enddo
ijkb0 = ijkb0 + nh (nt)
endif
enddo
else
do na = 1, nat
if (ityp (na) .eq.nt) ijkb0 = ijkb0 + nh (nt)
enddo
endif
enddo
!
deallocate (dbecq)
call stop_clock ('addusdbec')
return
end subroutine addusdbec