2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! 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 .
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine addusdbec (ik, wgt, psi, dbecsum)
|
2003-01-20 05:58:50 +08:00
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! 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.
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-06-12 21:44:18 +08:00
|
|
|
USE ions_base, ONLY : nat, ityp, ntyp => nsp
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2004-06-01 01:55:33 +08:00
|
|
|
USE uspp_param, only: nh, tvanp, nhm
|
2003-01-20 05:58:50 +08:00
|
|
|
use phcom
|
2003-02-08 00:04:36 +08:00
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! the dummy variables
|
|
|
|
!
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: dbecsum (nhm*(nhm+1)/2, nat), psi(npwx,nbnd)
|
2003-01-20 05:58:50 +08:00
|
|
|
! inp/out: the sum kv of bec *
|
|
|
|
! input : contains delta psi
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: ik
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the k point
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: wgt
|
2003-01-20 05:58:50 +08:00
|
|
|
! 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
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: w
|
2003-01-20 05:58:50 +08:00
|
|
|
! the weight of the k point
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: dbecq (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the change of becq
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
if (.not.okvan) return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call start_clock ('addusdbec')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (dbecq( nkb, nbndx))
|
|
|
|
w = 0.5d0 * wgt * omega
|
|
|
|
if (lgamma) then
|
|
|
|
ikk = ik
|
|
|
|
else
|
|
|
|
ikk = 2 * ik - 1
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
!
|
|
|
|
! First compute the product of psi and vkb
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call ccalbec (nkb, npwx, npwq, nbnd, dbecq, vkb, psi)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! And then we add the product to becsum
|
|
|
|
!
|
2003-10-15 19:40:07 +08:00
|
|
|
! Band parallelization: each processor takes care of its slice of bands
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call divide (nbnd_occ (ikk), startb, lastb)
|
2003-10-15 19:40:07 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
ijkb0 = 0
|
|
|
|
do nt = 1, ntyp
|
|
|
|
if (tvanp (nt) ) then
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) .eq.nt) then
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! And qgmq and becp and dbecq
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
ijh = 1
|
|
|
|
do ih = 1, nh (nt)
|
|
|
|
ikb = ijkb0 + ih
|
|
|
|
do ibnd = startb, lastb
|
2003-01-20 05:58:50 +08:00
|
|
|
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) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
ijh = ijh + 1
|
|
|
|
do jh = ih + 1, nh (nt)
|
|
|
|
jkb = ijkb0 + jh
|
|
|
|
do ibnd = startb, lastb
|
2003-01-20 05:58:50 +08:00
|
|
|
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) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
ijh = ijh + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
ijkb0 = ijkb0 + nh (nt)
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) .eq.nt) ijkb0 = ijkb0 + nh (nt)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
2003-10-15 19:40:07 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (dbecq)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call stop_clock ('addusdbec')
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine addusdbec
|