quantum-espresso/PW/force_ew.f90

170 lines
5.0 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 force_ew (alat, nat, ntyp, ityp, zv, at, bg, tau, &
omega, g, gg, ngm, gstart, gamma_only, gcutm, strf, forceion)
!-----------------------------------------------------------------------
!
! This routine computes the Ewald contribution to the forces,
! both the real- and reciprocal-space terms are present
!
#include "f_defs.h"
USE kinds
USE constants, ONLY : tpi, e2
USE mp_global, ONLY : intra_pool_comm
USE mp, ONLY : mp_sum
implicit none
!
! First the dummy variables
!
integer :: nat, ntyp, ngm, ityp (nat), gstart
! input: the number of atoms
! input: the number of types of atom
! input: the number of G vectors
! input: the type of each atom
! input: first non-zero G vector
logical :: gamma_only
real(DP) :: factor, tau (3, nat), g (3, ngm), gg (ngm), zv (ntyp), &
at (3, 3), bg (3, 3), omega, gcutm, alat
! input: the coordinates of the atoms
! input: the G vectors
! input: the moduli of G vectors
! input: the charge of the atoms
! input: the direct lattice vectors
! input: the reciprocal lattice vectors
! input: the volume of the unit cell
! input: cut-off of g vectors
! input: the edge of the cell
!
complex(DP) :: strf (ngm, ntyp)
! input: the structure factor on the potential
!
real(DP) :: forceion (3, nat)
! output: the ewald part of the forces
!
integer, parameter :: mxr=50
! the maximum number of R vectors
integer :: ig, n, na, nb, nt, nrm, ipol
! counter on G vectos
! counter on r vectors
! counter on atoms
! counter on atoms
! counter on atomic types
! the number of R vectors for real space su
! counter on polarization
real(DP) :: sumnb, arg, tpiba2, alpha, dtau (3), r (3, mxr), &
r2 (mxr), rmax, rr, charge, upperbound, fact
! auxiliary variable for speed
! the argument of the exponential
! 2 pi /alat
! the alpha parameter
! the difference of two tau
! the position of the atoms in the shell
! the square of r
! the maximum r
! the modulus of the r vectors
! the total charge
! used to determine alpha
complex(DP), allocatable :: aux (:)
! auxiliary space
real(DP), external :: erfc
!
forceion(:,:) = 0.d0
tpiba2 = (tpi / alat) **2
charge = 0.d0
do na = 1, nat
charge = charge+zv (ityp (na) )
enddo
!
! choose alpha in order to have convergence in the sum over G
! upperbound is a safe upper bound for the error ON THE ENERGY
!
alpha = 1.1d0
10 alpha = alpha - 0.1d0
if (alpha.eq.0.d0) call errore ('force_ew', 'optimal alpha not found', 1)
upperbound = e2 * charge**2 * sqrt (2.d0 * alpha / tpi) * &
erfc ( sqrt (tpiba2 * gcutm / 4.d0 / alpha) )
if (upperbound > 1.0d-6) goto 10
!
! G-space sum here
!
allocate(aux(ngm))
aux(:) = (0.d0, 0.d0)
do nt = 1, ntyp
do ig = 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
aux (ig) = aux (ig) + zv (nt) * CONJG(strf (ig, nt) )
enddo
enddo
do ig = gstart, ngm
aux (ig) = aux (ig) * exp ( - gg (ig) * tpiba2 / alpha / 4.d0) &
/ (gg (ig) * tpiba2)
enddo
if (gamma_only) then
fact = 4.d0
else
fact = 2.d0
end if
do na = 1, nat
do ig = gstart, ngm
arg = tpi * (g (1, ig) * tau (1, na) + g (2, ig) * tau (2, na) &
+ g (3, ig) * tau (3, na) )
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
sumnb = cos (arg) * AIMAG (aux(ig)) - sin (arg) * DBLE (aux(ig) )
forceion (1, na) = forceion (1, na) + g (1, ig) * sumnb
forceion (2, na) = forceion (2, na) + g (2, ig) * sumnb
forceion (3, na) = forceion (3, na) + g (3, ig) * sumnb
enddo
do ipol = 1, 3
forceion (ipol, na) = - zv (ityp (na) ) * fact * e2 * tpi**2 / &
omega / alat * forceion (ipol, na)
enddo
enddo
deallocate (aux)
if (gstart == 1) goto 100
!
! R-space sum here (only for the processor that contains G=0)
!
rmax = 5.d0 / (sqrt (alpha) * alat)
!
! with this choice terms up to ZiZj*erfc(5) are counted (erfc(5)=2x10^-1
!
do na = 1, nat
do nb = 1, nat
if (nb.eq.na) goto 50
dtau (:) = tau (:, na) - tau (:, nb)
!
! generates nearest-neighbors shells r(i)=R(i)-dtau(i)
!
call rgen (dtau, rmax, mxr, at, bg, r, r2, nrm)
do n = 1, nrm
rr = sqrt (r2 (n) ) * alat
factor = zv (ityp (na) ) * zv (ityp (nb) ) * e2 / rr**2 * &
(erfc (sqrt (alpha) * rr) / rr + &
sqrt (8.0d0 * alpha / tpi) * exp ( - alpha * rr**2) ) * alat
do ipol = 1, 3
forceion (ipol, na) = forceion (ipol, na) - factor * r (ipol, n)
enddo
enddo
50 continue
enddo
enddo
100 continue
#ifdef __PARA
call mp_sum( forceion, intra_pool_comm )
#endif
return
end subroutine force_ew