2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-03-14 02:20:45 +08:00
|
|
|
! Copyright (C) 2001-2003 PWSCF 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 .
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine set_rhoc
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
2003-03-14 02:20:45 +08:00
|
|
|
! This routine computes the core charge on the real space 3D mesh
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
!
|
2004-06-12 21:44:18 +08:00
|
|
|
USE io_global, ONLY : stdout
|
|
|
|
USE kinds, ONLY : DP
|
|
|
|
USE atom, ONLY : rho_atc, numeric, msh, r, rab, nlcc
|
|
|
|
USE ions_base, ONLY : ntyp => nsp
|
|
|
|
USE cell_base, ONLY : omega, tpiba2
|
|
|
|
USE ener, ONLY : etxcc
|
|
|
|
USE gvect, ONLY : ngm, nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
|
|
|
nrxx, nl, nlm, ngl, gl, igtongl
|
|
|
|
USE pseud, ONLY : a_nlcc, b_nlcc, alpha_nlcc
|
|
|
|
USE scf, ONLY : rho_core
|
|
|
|
USE vlocal, ONLY : strf
|
|
|
|
USE wvfct, ONLY : gamma_only
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
implicit none
|
|
|
|
!
|
2003-03-14 02:20:45 +08:00
|
|
|
real(kind=DP), parameter :: eps = 1.d-10
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP) , allocatable :: aux (:)
|
2003-03-14 02:20:45 +08:00
|
|
|
! used for the fft of the core charge
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
real(kind=DP) , allocatable :: rhocg(:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the radial fourier trasform
|
2003-03-14 02:20:45 +08:00
|
|
|
real(kind=DP) :: rhoima, rhoneg, rhorea
|
2003-01-20 05:58:50 +08:00
|
|
|
! used to check the core charge
|
2003-03-14 02:20:45 +08:00
|
|
|
real(kind=DP) :: vtxcc
|
2003-01-20 05:58:50 +08:00
|
|
|
! dummy xc energy term
|
2003-03-14 02:20:45 +08:00
|
|
|
real(kind=DP) , allocatable :: dum(:,:)
|
|
|
|
! dummy array containing rho=0
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: ir, nt, ng
|
2003-01-20 05:58:50 +08:00
|
|
|
! counter on mesh points
|
|
|
|
! counter on atomic types
|
|
|
|
! counter on g vectors
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
etxcc = 0.d0
|
|
|
|
do nt = 1, ntyp
|
|
|
|
if (nlcc (nt) ) goto 10
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-04-23 00:03:45 +08:00
|
|
|
rho_core(:) = 0.d0
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
10 continue
|
|
|
|
allocate (aux( nrxx))
|
2003-04-23 00:03:45 +08:00
|
|
|
allocate (rhocg( ngl))
|
2003-03-14 02:20:45 +08:00
|
|
|
aux (:) = 0.d0
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! the sum is on atom types
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do nt = 1, ntyp
|
|
|
|
if (nlcc (nt) ) then
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! drhoc compute the radial fourier transform for each shell of g vec
|
|
|
|
!
|
|
|
|
call drhoc (ngl, gl, omega, tpiba2, numeric (nt), a_nlcc (nt), &
|
|
|
|
b_nlcc (nt), alpha_nlcc (nt), msh (nt), r (1, nt), rab (1, nt), &
|
|
|
|
rho_atc (1, nt), rhocg)
|
|
|
|
!
|
|
|
|
! multiply by the structure factor and sum
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do ng = 1, ngm
|
2003-04-23 00:03:45 +08:00
|
|
|
aux(nl(ng)) = aux(nl(ng)) + strf(ng,nt) * rhocg(igtongl(ng))
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
2003-10-30 02:53:40 +08:00
|
|
|
if (gamma_only) then
|
|
|
|
do ng = 1, 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
|
|
|
aux(nlm(ng)) = CONJG(aux(nl (ng)))
|
2003-10-30 02:53:40 +08:00
|
|
|
end do
|
|
|
|
end if
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! the core charge in real space
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call cft3 (aux, nr1, nr2, nr3, nrx1, nrx2, nrx3, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! test on the charge and computation of the core energy
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
rhoneg = 0.d0
|
|
|
|
rhoima = 0.d0
|
|
|
|
do ir = 1, nrxx
|
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
|
|
|
rhoneg = rhoneg + min (0.d0, DBLE (aux (ir) ) )
|
|
|
|
rhoima = rhoima + abs (AIMAG (aux (ir) ) )
|
|
|
|
rho_core(ir) = DBLE (aux(ir))
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! NOTE: Core charge is computed in reciprocal space and brought to real
|
|
|
|
! space by FFT. For non smooth core charges (or insufficient cut-off)
|
|
|
|
! this may result in negative values in some grid points.
|
|
|
|
! Up to October 1999 the core charge was forced to be positive definite.
|
2003-04-23 00:03:45 +08:00
|
|
|
! This induces an error in the force, and probably stress, calculation if
|
2003-01-20 05:58:50 +08:00
|
|
|
! the number of grid points where the core charge would be otherwise neg
|
2003-04-23 00:03:45 +08:00
|
|
|
! is large. The error disappears for sufficiently high cut-off, but may be
|
2003-01-20 05:58:50 +08:00
|
|
|
! rather large and it is better to leave the core charge as it is.
|
|
|
|
! If you insist to have it positive definite (with the possible problems
|
2003-03-14 02:20:45 +08:00
|
|
|
! mentioned above) uncomment the following lines. SdG, Oct 15 1999
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
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
|
|
|
! rhorea = max ( DBLE (aux (ir) ), eps)
|
2003-01-20 05:58:50 +08:00
|
|
|
! rho_core(ir) = rhorea
|
|
|
|
!
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
rhoneg = rhoneg / (nr1 * nr2 * nr3)
|
|
|
|
rhoima = rhoima / (nr1 * nr2 * nr3)
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2003-02-08 00:04:36 +08:00
|
|
|
call reduce (1, rhoneg)
|
|
|
|
call reduce (1, rhoima)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
2005-03-05 01:49:28 +08:00
|
|
|
IF (rhoneg < -1.0d-6 .OR. rhoima > 1.0d-6) &
|
|
|
|
WRITE( stdout, '(/5x,"warning: negative or imaginary core charge ",2f12.6)')&
|
2003-01-20 05:58:50 +08:00
|
|
|
rhoneg, rhoima
|
|
|
|
!
|
2003-03-14 02:20:45 +08:00
|
|
|
! calculate core_only exch-corr energy etxcc=E_xc[rho_core] if required
|
|
|
|
! The term was present in previous versions of the code but it shouldn't
|
|
|
|
!
|
|
|
|
! allocate (dum(nrxx , nspin))
|
|
|
|
! dum(:,:) = 0.d0
|
|
|
|
! call v_xc (dum, rho_core, nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
2003-12-17 23:38:44 +08:00
|
|
|
! nrxx, nl, ngm, g, nspin, alat, omega, etxcc, vtxcc, aux)
|
|
|
|
!
|
2003-03-14 02:20:45 +08:00
|
|
|
! deallocate(dum)
|
2003-11-04 18:53:05 +08:00
|
|
|
! WRITE( stdout, 9000) etxcc
|
|
|
|
! WRITE( stdout, * ) 'BEWARE it will be subtracted from total energy !'
|
2003-03-14 02:20:45 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (rhocg)
|
|
|
|
deallocate (aux)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
9000 format (5x,'core-only xc energy = ',f15.8,' ryd')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
end subroutine set_rhoc
|
|
|
|
|