2003-01-20 05:58:50 +08:00
|
|
|
!
|
2007-12-19 06:31:46 +08:00
|
|
|
! Copyright (C) 2001-2007 Quantum-ESPRESSO group
|
2003-01-20 05:58:50 +08:00
|
|
|
! 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 drhodv (nu_i0, nper, drhoscf)
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This subroutine computes the electronic term
|
|
|
|
! <psi|dv|dpsi> of the dynamical matrix
|
|
|
|
!
|
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
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
2007-12-19 06:31:46 +08:00
|
|
|
USE becmod, ONLY: calbec
|
2007-02-08 21:40:53 +08:00
|
|
|
USE io_global, ONLY : stdout
|
|
|
|
USE noncollin_module, ONLY : noncolin, npol
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2003-11-10 02:30:08 +08:00
|
|
|
USE io_files, ONLY: iunigk
|
2003-02-08 00:04:36 +08:00
|
|
|
use phcom
|
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: nper, nu_i0
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: number of perturbations of this represent
|
|
|
|
! input: the initial position of the mode
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: drhoscf (nrxx, nspin, npertx)
|
2004-03-07 21:47:42 +08:00
|
|
|
! the change of density due to perturbations
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
integer :: mu, ik, ikq, ig, nu_i, nu_j, na_jcart, ibnd, nrec, &
|
|
|
|
ipol, ikk
|
2004-03-07 21:47:42 +08:00
|
|
|
! counters
|
|
|
|
! ikk: record position for wfc at k
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: fact, ps, dynwrk (3 * nat, 3 * nat), &
|
2004-03-07 21:47:42 +08:00
|
|
|
wdyn (3 * nat, 3 * nat), ZDOTC
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: aux (:,:), dbecq (:,:,:), &
|
2007-02-08 21:40:53 +08:00
|
|
|
dalpq (:,:,:,:), dbecq_nc(:,:,:,:), dalpq_nc(:,:,:,:,:)
|
2004-03-07 21:47:42 +08:00
|
|
|
! work space
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Initialize the auxiliary matrix wdyn
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call start_clock ('drhodv')
|
2007-02-08 21:40:53 +08:00
|
|
|
if (noncolin) then
|
|
|
|
allocate (dbecq_nc ( nkb, npol, nbnd, nper))
|
|
|
|
allocate (dalpq_nc ( nkb, npol, nbnd ,3 ,nper))
|
|
|
|
else
|
|
|
|
allocate (dbecq ( nkb , nbnd, nper))
|
|
|
|
allocate (dalpq ( nkb , nbnd ,3 ,nper))
|
|
|
|
endif
|
|
|
|
allocate (aux ( npwx*npol , nbnd))
|
2004-03-07 21:47:42 +08:00
|
|
|
dynwrk(:,:) = (0.d0, 0.d0)
|
|
|
|
wdyn (:,:) = (0.d0, 0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! We need a sum over all k points ...
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nksq > 1) rewind (unit = iunigk)
|
2003-02-08 00:04:36 +08:00
|
|
|
do ik = 1, nksq
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nksq > 1) read (iunigk) npw, igk
|
2003-02-08 00:04:36 +08:00
|
|
|
if (lgamma) then
|
|
|
|
ikk = ik
|
|
|
|
ikq = ik
|
|
|
|
npwq = npw
|
|
|
|
else
|
|
|
|
ikk = 2 * ik - 1
|
|
|
|
ikq = ikk + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nksq > 1) read (iunigk) npwq, igkq
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
if (lsda) current_spin = isk (ikk)
|
|
|
|
call init_us_2 (npwq, igkq, xk (1, ikq), vkb)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
do mu = 1, nper
|
|
|
|
nrec = (mu - 1) * nksq + ik
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nksq > 1 .or. nper > 1) call davcio(dpsi, lrdwf, iudwf, nrec,-1)
|
2007-02-08 21:40:53 +08:00
|
|
|
if (noncolin) then
|
2007-12-19 06:31:46 +08:00
|
|
|
call calbec (npwq, vkb, dpsi, dbecq_nc(:,:,:,mu) )
|
2007-02-08 21:40:53 +08:00
|
|
|
else
|
2007-12-19 06:31:46 +08:00
|
|
|
call calbec (npwq, vkb, dpsi, dbecq(:,:,mu) )
|
2007-02-08 21:40:53 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipol = 1, 3
|
2007-02-08 21:40:53 +08:00
|
|
|
aux=(0.d0,0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
do ibnd = 1, nbnd
|
|
|
|
do ig = 1, npwq
|
2003-01-20 05:58:50 +08:00
|
|
|
aux (ig, ibnd) = dpsi (ig, ibnd) * &
|
|
|
|
(xk (ipol, ikq) + g (ipol, igkq (ig) ) )
|
|
|
|
enddo
|
2007-02-08 21:40:53 +08:00
|
|
|
if (noncolin) then
|
|
|
|
do ig = 1, npwq
|
|
|
|
aux (ig+npwx, ibnd) = dpsi (ig+npwx, ibnd) * &
|
|
|
|
(xk (ipol, ikq) + g (ipol, igkq (ig) ) )
|
|
|
|
enddo
|
|
|
|
endif
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2007-02-08 21:40:53 +08:00
|
|
|
if (noncolin) then
|
2007-12-19 06:31:46 +08:00
|
|
|
call calbec (npwq, vkb, aux, dalpq_nc(:,:,:,ipol,mu) )
|
2007-02-08 21:40:53 +08:00
|
|
|
else
|
2007-12-19 06:31:46 +08:00
|
|
|
call calbec (npwq, vkb, aux, dalpq(:,:,ipol,mu) )
|
2007-02-08 21:40:53 +08:00
|
|
|
endif
|
2003-01-20 05:58:50 +08:00
|
|
|
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)
|
2007-02-08 21:40:53 +08:00
|
|
|
if (noncolin) then
|
|
|
|
dalpq_nc = dalpq_nc * fact
|
|
|
|
call drhodvnl (ik, ikk, nper, nu_i0, dynwrk, dbecq_nc, dalpq_nc)
|
|
|
|
else
|
|
|
|
dalpq = dalpq * fact
|
|
|
|
call drhodvnl (ik, ikk, nper, nu_i0, dynwrk, dbecq, dalpq)
|
|
|
|
endif
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! put in the basis of the modes
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do nu_i = 1, 3 * nat
|
|
|
|
do nu_j = 1, 3 * nat
|
|
|
|
ps = (0.0d0, 0.0d0)
|
|
|
|
do na_jcart = 1, 3 * nat
|
|
|
|
ps = ps + dynwrk (nu_i, na_jcart) * u (na_jcart, nu_j)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
wdyn (nu_i, nu_j) = wdyn (nu_i, nu_j) + ps
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! collect contributions from all pools (sum over k-points)
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call poolreduce (18 * nat * nat, wdyn)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! add the contribution of the local part of the perturbation
|
|
|
|
!
|
2007-02-08 21:40:53 +08:00
|
|
|
call drhodvloc (nu_i0, nper, drhoscf, wdyn)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! add to the rest of the dynamical matrix
|
|
|
|
!
|
2003-11-06 03:01:20 +08:00
|
|
|
! WRITE( stdout,*) 'drhodv dyn, wdyn'
|
2003-01-20 05:58:50 +08:00
|
|
|
! call tra_write_matrix('drhodv dyn',dyn,u,nat)
|
|
|
|
! call tra_write_matrix('drhodv wdyn',wdyn,u,nat)
|
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
dyn (:,:) = dyn (:,:) + wdyn (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (aux)
|
2007-02-08 21:40:53 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
deallocate (dalpq_nc)
|
|
|
|
deallocate (dbecq_nc)
|
|
|
|
ELSE
|
|
|
|
deallocate (dalpq)
|
|
|
|
deallocate (dbecq)
|
|
|
|
ENDIF
|
2003-02-08 00:04:36 +08:00
|
|
|
|
|
|
|
call stop_clock ('drhodv')
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine drhodv
|