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 addnlcc (imode0, drhoscf, npe)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! This routine adds a contribution to the dynamical matrix due
|
|
|
|
! to the NLCC
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-06-12 21:44:18 +08:00
|
|
|
USE ions_base, ONLY : nat
|
2003-01-20 05:58:50 +08:00
|
|
|
use funct
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2003-01-20 05:58:50 +08:00
|
|
|
use phcom
|
|
|
|
implicit none
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: imode0, npe
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the starting mode
|
|
|
|
! input: the number of perturbations
|
|
|
|
! input: the change of density due to perturbatio
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP) :: drhoscf (nrxx, nspin, npertx)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: nrtot, ipert, jpert, is, is1, irr, ir, mode, mode1
|
2003-01-20 05:58:50 +08:00
|
|
|
! the total number of points
|
|
|
|
! counter on perturbations
|
|
|
|
! counter on spin
|
|
|
|
! counter on representations
|
|
|
|
! counter on real space points
|
|
|
|
! counter on modes
|
|
|
|
|
|
|
|
complex(kind=DP) :: ZDOTC, dyn1 (3 * nat, 3 * nat)
|
|
|
|
! the scalar product function
|
|
|
|
! auxiliary dynamical matrix
|
2004-03-07 21:47:42 +08:00
|
|
|
complex(kind=DP), allocatable :: drhoc (:), dvaux (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the change of the core
|
|
|
|
! the change of the potential
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
real(kind=DP) :: fac
|
2003-01-20 05:58:50 +08:00
|
|
|
! auxiliary factor
|
|
|
|
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
if (.not.nlcc_any) return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (drhoc( nrxx))
|
|
|
|
allocate (dvaux( nrxx , nspin))
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
dyn1 (:,:) = (0.d0, 0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! compute the exchange and correlation potential for this mode
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
nrtot = nr1 * nr2 * nr3
|
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
|
|
|
fac = 1.d0 / DBLE (nspin)
|
2003-02-08 00:04:36 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
do ipert = 1, npe
|
2003-02-08 00:04:36 +08:00
|
|
|
mode = imode0 + ipert
|
2004-03-07 21:47:42 +08:00
|
|
|
dvaux (:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
call addcore (mode, drhoc)
|
|
|
|
do is = 1, nspin
|
|
|
|
call DAXPY (nrxx, fac, rho_core, 1, rho (1, is), 1)
|
2003-10-16 22:39:25 +08:00
|
|
|
call DAXPY (2 * nrxx, fac, drhoc, 1, drhoscf (1, is, ipert), 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
|
|
|
do is1 = 1, nspin
|
|
|
|
do ir = 1, nrxx
|
2003-10-16 22:39:25 +08:00
|
|
|
dvaux (ir, is) = dvaux (ir, is) + dmuxc (ir, is, is1) * &
|
|
|
|
drhoscf ( ir, is1, ipert)
|
2003-01-20 05:58:50 +08:00
|
|
|
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
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
if (igcx /= 0 .or. igcc /= 0) &
|
2003-01-20 05:58:50 +08:00
|
|
|
call dgradcorr (rho, grho, dvxc_rr, dvxc_sr, dvxc_ss, dvxc_s, xq, &
|
|
|
|
drhoscf (1, 1, ipert), nr1, nr2, nr3, nrx1, nrx2, nrx3, nrxx, &
|
|
|
|
nspin, nl, ngm, g, alat, omega, dvaux)
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
|
|
|
call DAXPY (nrxx, - fac, rho_core, 1, rho (1, is), 1)
|
2003-10-16 22:39:25 +08:00
|
|
|
call DAXPY (2 * nrxx, - fac, drhoc, 1, drhoscf (1, is, ipert), 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
mode1 = 0
|
|
|
|
do irr = 1, nirr
|
|
|
|
do jpert = 1, npert (irr)
|
|
|
|
mode1 = mode1 + 1
|
|
|
|
call addcore (mode1, drhoc)
|
|
|
|
do is = 1, nspin
|
2004-03-07 21:47:42 +08:00
|
|
|
dyn1 (mode, mode1) = dyn1 (mode, mode1) + &
|
|
|
|
ZDOTC (nrxx, dvaux (1, is), 1, drhoc, 1) * &
|
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
|
|
|
omega * fac / DBLE (nrtot)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
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 reduce (18 * nat * nat, dyn1)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
2004-03-07 21:47:42 +08:00
|
|
|
dyn (:,:) = dyn(:,:) + dyn1(:,:)
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (dvaux)
|
|
|
|
deallocate (drhoc)
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine addnlcc
|