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 dv_of_drho (mode, dvscf, flag)
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine computes the change of the self consistent potential
|
|
|
|
! due to the perturbation.
|
|
|
|
!
|
|
|
|
#include "machine.h"
|
|
|
|
use funct
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
|
|
|
use parameters, only : DP
|
|
|
|
use phcom
|
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: mode
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the mode to do
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP) :: dvscf (nrxx, nspin)
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the change of the charge,
|
|
|
|
! output: change of the potential
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
logical :: flag
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: if true add core charge
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: ir, is, is1, ig
|
2003-01-20 05:58:50 +08:00
|
|
|
! counter on r vectors
|
|
|
|
! counter on spin polarizations
|
|
|
|
! counter on g vectors
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
real(kind=DP) :: qg2, fac
|
2003-01-20 05:58:50 +08:00
|
|
|
! the modulus of (q+G)^2
|
|
|
|
! the structure factor
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP), allocatable :: dvaux (:,:), drhoc (:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! auxiliary variable for potential
|
|
|
|
! the change of the core charge
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call start_clock ('dv_of_drho')
|
|
|
|
allocate (dvaux( nrxx, nspin))
|
|
|
|
allocate (drhoc( nrxx))
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! the exchange-correlation contribution is computed in real space
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call setv (2 * nrxx * nspin, 0.d0, dvaux, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
fac = 1.d0 / float (nspin)
|
|
|
|
if (nlcc_any.and.flag) then
|
|
|
|
call addcore (mode, drhoc)
|
|
|
|
do is = 1, nspin
|
|
|
|
call DAXPY (nrxx, fac, rho_core, 1, rho (1, is), 1)
|
|
|
|
call DAXPY (2 * nrxx, fac, drhoc, 1, dvscf (1, is), 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
|
|
|
do is1 = 1, nspin
|
|
|
|
do ir = 1, nrxx
|
2003-01-20 05:58:50 +08:00
|
|
|
dvaux (ir, is) = dvaux (ir, is) + dmuxc (ir, is, is1) * dvscf (ir, &
|
|
|
|
is1)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! add gradient correction to xc, NB: if nlcc is true we need to add here
|
|
|
|
! its contribution. grho contains already the core charge
|
|
|
|
!
|
|
|
|
if (igcx.ne.0.or.igcc.ne.0) call dgradcorr &
|
|
|
|
(rho, grho, dvxc_rr, dvxc_sr, dvxc_ss, dvxc_s, xq, &
|
|
|
|
dvscf, nr1, nr2, nr3, nrx1, nrx2, nrx3, nrxx, nspin, nl, ngm, g, &
|
|
|
|
alat, omega, dvaux)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (nlcc_any.and.flag) then
|
|
|
|
do is = 1, nspin
|
|
|
|
call DAXPY (nrxx, - fac, rho_core, 1, rho (1, is), 1)
|
|
|
|
call DAXPY (2 * nrxx, - fac, drhoc, 1, dvscf (1, is), 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
! copy the total (up+down) delta rho in dvscf(*,1) and go to G-space
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 2, nspin
|
2003-01-20 05:58:50 +08:00
|
|
|
call DAXPY (2 * nrxx, 1.d0, dvscf (1, is), 1, dvscf (1, 1), &
|
|
|
|
1)
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
call cft3 (dvscf, nr1, nr2, nr3, nrx1, nrx2, nrx3, - 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! hartree contribution is computed in reciprocal space
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
2003-01-20 05:58:50 +08:00
|
|
|
call cft3 (dvaux (1, is), nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
|
|
|
- 1)
|
2003-02-08 00:04:36 +08:00
|
|
|
do ig = 1, ngm
|
2003-01-20 05:58:50 +08:00
|
|
|
qg2 = (g (1, ig) + xq (1) ) **2 + (g (2, ig) + xq (2) ) **2 + &
|
|
|
|
(g (3, ig) + xq (3) ) **2
|
2003-02-08 00:04:36 +08:00
|
|
|
if (qg2.gt.1.d-8) then
|
2003-01-20 05:58:50 +08:00
|
|
|
dvaux (nl (ig), is) = dvaux (nl (ig), is) + e2 * fpi * dvscf ( &
|
|
|
|
nl (ig), 1) / (tpiba2 * qg2)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! and transformed back to real space
|
|
|
|
!
|
|
|
|
call cft3 (dvaux (1, is), nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
|
|
|
+ 1)
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! at the end the two contributes are added
|
|
|
|
!
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call ZCOPY (nrxx * nspin, dvaux, 1, dvscf, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (drhoc)
|
|
|
|
deallocate (dvaux)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call stop_clock ('dv_of_drho')
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine dv_of_drho
|