quantum-espresso/PH/addnlcc.f90

112 lines
3.1 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 addnlcc (imode0, drhoscf, npe)
!
! This routine adds a contribution to the dynamical matrix due
! to the NLCC
!
#include "f_defs.h"
USE ions_base, ONLY : nat
use funct
use pwcom
USE kinds, only : DP
use phcom
implicit none
integer :: imode0, npe
! input: the starting mode
! input: the number of perturbations
! input: the change of density due to perturbatio
complex(kind=DP) :: drhoscf (nrxx, nspin, npertx)
integer :: nrtot, ipert, jpert, is, is1, irr, ir, mode, mode1
! 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
complex(kind=DP), allocatable :: drhoc (:), dvaux (:,:)
! the change of the core
! the change of the potential
real(kind=DP) :: fac
! auxiliary factor
if (.not.nlcc_any) return
allocate (drhoc( nrxx))
allocate (dvaux( nrxx , nspin))
dyn1 (:,:) = (0.d0, 0.d0)
!
! compute the exchange and correlation potential for this mode
!
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)
do ipert = 1, npe
mode = imode0 + ipert
dvaux (:,:) = (0.d0, 0.d0)
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, drhoscf (1, is, ipert), 1)
enddo
do is = 1, nspin
do is1 = 1, nspin
do ir = 1, nrxx
dvaux (ir, is) = dvaux (ir, is) + dmuxc (ir, is, is1) * &
drhoscf ( ir, is1, ipert)
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 /= 0 .or. igcc /= 0) &
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)
do is = 1, nspin
call DAXPY (nrxx, - fac, rho_core, 1, rho (1, is), 1)
call DAXPY (2 * nrxx, - fac, drhoc, 1, drhoscf (1, is, ipert), 1)
enddo
mode1 = 0
do irr = 1, nirr
do jpert = 1, npert (irr)
mode1 = mode1 + 1
call addcore (mode1, drhoc)
do is = 1, nspin
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)
enddo
enddo
enddo
enddo
#ifdef __PARA
!
! collect contributions from all pools (sum over k-points)
!
call reduce (18 * nat * nat, dyn1)
#endif
dyn (:,:) = dyn(:,:) + dyn1(:,:)
deallocate (dvaux)
deallocate (drhoc)
return
end subroutine addnlcc