quantum-espresso/PW/stres_cc.f90

112 lines
3.3 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 .
!
#include "f_defs.h"
!
!-----------------------------------------------------------------------
subroutine stres_cc (sigmaxcc)
!-----------------------------------------------------------------------
!
USE kinds, ONLY : DP
USE atom, ONLY : rho_atc, numeric, mesh, r, rab, nlcc
USE ions_base, ONLY : ntyp => nsp
USE cell_base, ONLY : alat, omega, tpiba, tpiba2
USE gvect, ONLY : ngm, gstart, nr1, nr2, nr3, nrx1, nrx2, &
nrx3, nrxx, nl, g, gg, ngl, gl, igtongl
USE ener, ONLY : etxc, vtxc
USE lsda_mod, ONLY : nspin
USE pseud, ONLY : a_nlcc, b_nlcc, alpha_nlcc
USE scf, ONLY : rho, rho_core
USE vlocal, ONLY : strf
USE wvfct, ONLY : gamma_only
USE wavefunctions_module, ONLY : psic
!
implicit none
! output
real(DP) :: sigmaxcc (3, 3)
! local variables
integer :: nt, ng, l, m, ir
! counters
real(DP) :: fact, sigmadiag
real(DP) , allocatable:: rhocg (:), vxc (:,:)
sigmaxcc(:,:) = 0.d0
do nt = 1, ntyp
if (nlcc (nt) ) goto 15
enddo
return
15 continue
!
! recalculate the exchange-correlation potential
!
allocate ( vxc(nrxx,nspin) )
call v_xc (rho, rho_core, nr1, nr2, nr3, nrx1, nrx2, nrx3, nrxx, &
nl, ngm, g, nspin, alat, omega, etxc, vtxc, vxc)
if (nspin.eq.1) then
do ir = 1, nrxx
psic (ir) = vxc (ir, 1)
enddo
else
do ir = 1, nrxx
psic (ir) = 0.5d0 * (vxc (ir, 1) + vxc (ir, 2) )
enddo
endif
deallocate (vxc)
call cft3 (psic, nr1, nr2, nr3, nrx1, nrx2, nrx3, - 1)
!
! psic contains now Vxc(G)
!
allocate(rhocg(ngl))
sigmadiag = 0.0
if (gamma_only) then
fact = 2.d0
else
fact = 1.d0
end if
do nt = 1, ntyp
if (nlcc (nt) ) then
call drhoc (ngl, gl, omega, tpiba2, numeric (nt), a_nlcc (nt), &
b_nlcc (nt), alpha_nlcc (nt), mesh (nt), r (1, nt), rab (1, nt) &
, rho_atc (1, nt), rhocg)
! diagonal term
if (gstart==2) sigmadiag = sigmadiag + &
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
CONJG(psic (nl(1) ) ) * strf (1,nt) * rhocg (igtongl (1) )
do ng = gstart, ngm
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
sigmadiag = sigmadiag + CONJG(psic (nl (ng) ) ) * &
strf (ng,nt) * rhocg (igtongl (ng) ) * fact
enddo
call deriv_drhoc (ngl, gl, omega, tpiba2, numeric (nt), &
a_nlcc (nt), b_nlcc (nt), alpha_nlcc (nt), mesh (nt), r (1, nt), &
rab (1, nt), rho_atc (1, nt), rhocg)
! non diagonal term (g=0 contribution missing)
do ng = gstart, ngm
do l = 1, 3
do m = 1, 3
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
sigmaxcc (l, m) = sigmaxcc (l, m) + CONJG(psic (nl (ng) ) ) &
* strf (ng, nt) * rhocg (igtongl (ng) ) * tpiba * &
g (l, ng) * g (m, ng) / sqrt (gg (ng) ) * fact
enddo
enddo
enddo
endif
enddo
do l = 1, 3
sigmaxcc (l, l) = sigmaxcc (l, l) + sigmadiag
enddo
#ifdef __PARA
call reduce (9, sigmaxcc)
#endif
deallocate (rhocg)
return
end subroutine stres_cc