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 .
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
|
|
|
|
subroutine d2ionq (nat, ntyp, ityp, zv, tau, alat, omega, q, at, &
|
|
|
|
bg, g, gg, ngm, gcutm, nmodes, u, dyn)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine computes the contribution of the ions to the
|
|
|
|
! dynamical matrix. Both the real and reciprocal space terms
|
|
|
|
! are included.
|
|
|
|
!
|
|
|
|
! The original routine was from C. Bungaro.
|
|
|
|
! Revised 16 oct. 1995 by Andrea Dal Corso.
|
|
|
|
! April 1997: parallel stuff added (SdG)
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2003-11-06 03:01:20 +08:00
|
|
|
USE io_global, ONLY : stdout
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2004-03-07 21:47:42 +08:00
|
|
|
USE constants, ONLY: e2, tpi, fpi
|
2008-04-19 01:05:56 +08:00
|
|
|
USE mp_global, ONLY: intra_pool_comm
|
|
|
|
USE mp, ONLY: mp_sum
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! first the dummy variables
|
|
|
|
!
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: nat, ntyp, ngm, ityp (nat), nmodes
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the number of atoms
|
|
|
|
! input: the number of types of atoms
|
|
|
|
! input: the number of G vectors
|
|
|
|
! input: the type of each atom
|
|
|
|
! input: the number of modes
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: tau (3, nat), g (3, ngm), gg (ngm), zv (ntyp), &
|
2003-01-20 05:58:50 +08:00
|
|
|
at (3, 3), bg (3, 3), omega, alat, gcutm, q (3)
|
|
|
|
! input: the positions of the atoms
|
|
|
|
! input: the coordinates of g vectors
|
|
|
|
! input: the modulus of g vectors
|
|
|
|
! input: the charge of each type
|
|
|
|
! input: the direct lattice vectors
|
|
|
|
! input: the reciprocal lattice vectors
|
|
|
|
! input: the volume of the unit cell
|
|
|
|
! input: the length scale
|
|
|
|
! input: cut-off of g vectors
|
|
|
|
! input: the q vector
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: dyn (3 * nat, nmodes), u (3 * nat, nmodes)
|
2003-01-20 05:58:50 +08:00
|
|
|
! output: the ionic part of the dyn. mat
|
|
|
|
! input: the pattern of the modes
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
! Local variables
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
integer, parameter :: mxr = 100
|
2003-01-20 05:58:50 +08:00
|
|
|
! the maximum number of r shells
|
|
|
|
|
|
|
|
integer :: nu_i, nu_j, na, nb, nta, ntb, ng, nrm, nr, icart, &
|
|
|
|
jcart, na_icart, na_jcart, nb_icart, nb_jcart
|
2004-03-07 21:47:42 +08:00
|
|
|
! counters
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: arg, argq, tpiba2, alpha, r (3, mxr), r2 (mxr), &
|
2009-02-24 01:44:13 +08:00
|
|
|
dtau (3), rmax, rr, upperbound, charge, fac, df, d2f, ar, &
|
2003-01-20 05:58:50 +08:00
|
|
|
gtq2, gt2, facq, qrg
|
2004-03-07 21:47:42 +08:00
|
|
|
! auxiliary variables
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: dy1 (3 * nat, nmodes), dy2 (3 * nat, nmodes), &
|
2003-01-20 05:58:50 +08:00
|
|
|
dy3 (3 * nat, nmodes), facg, fnat, work
|
2004-03-07 21:47:42 +08:00
|
|
|
! work spaces, factors
|
2009-02-24 01:44:13 +08:00
|
|
|
real(DP), external :: erfc
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call start_clock ('d2ionq')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
tpiba2 = (tpi / alat) **2
|
|
|
|
charge = 0.d0
|
|
|
|
do na = 1, nat
|
|
|
|
charge = charge+zv (ityp (na) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! choose alpha in order to have convergence in the sum over G
|
|
|
|
! upperbound is an upper bound for the error in the sum over G
|
|
|
|
! estimated for the energy (empirical trust!)
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
alpha = 2.9d0
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
11 alpha = alpha - 0.1d0
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (alpha == 0.d0) call errore ('d2ion', 'optimal alpha not found',1)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
upperbound = 2.d0 * charge**2 * sqrt (2.d0 * alpha / tpi) * &
|
|
|
|
erfc ( sqrt (tpiba2 * gcutm / 4.d0 / alpha) )
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (upperbound > 1.d-9) goto 11
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-11-06 03:01:20 +08:00
|
|
|
WRITE( stdout, '(/5x,"Alpha used in Ewald sum = ",f8.4)') alpha
|
2004-03-07 21:47:42 +08:00
|
|
|
dy1 (:,:) = (0.d0, 0.d0)
|
|
|
|
dy2 (:,:) = (0.d0, 0.d0)
|
|
|
|
dy3 (:,:) = (0.d0, 0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! G-space sums here
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do ng = 1, ngm
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! The reciprocal space term has two contributions. The first one
|
|
|
|
!
|
|
|
|
gtq2 = ( (g (1, ng) + q (1) ) **2 + (g (2, ng) + q (2) ) **2 + &
|
2004-03-07 21:47:42 +08:00
|
|
|
(g (3, ng) + q (3) ) **2) * tpiba2
|
|
|
|
if (abs (gtq2) > 1.d-8) then
|
|
|
|
facq = - e2*fpi * tpiba2 / omega * exp ( - gtq2 / alpha / 4.d0) / gtq2
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
|
|
|
facq = 0.d0
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
nta = ityp (na)
|
|
|
|
do nb = 1, nat
|
|
|
|
ntb = ityp (nb)
|
2003-01-20 05:58:50 +08:00
|
|
|
argq = tpi * ( (g (1, ng) + q (1) ) * (tau (1, na) - tau (1, nb) ) &
|
2004-03-07 21:47:42 +08:00
|
|
|
+ (g (2, ng) + q (2) ) * (tau (2, na) - tau (2, nb) ) &
|
|
|
|
+ (g (3, ng) + q (3) ) * (tau (3, na) - tau (3, nb) ) )
|
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
|
|
|
facg = facq * zv (nta) * zv (ntb) * CMPLX (cos (argq), sin (argq) )
|
2003-02-08 00:04:36 +08:00
|
|
|
do icart = 1, 3
|
|
|
|
nu_i = 3 * (na - 1) + icart
|
|
|
|
do jcart = 1, 3
|
|
|
|
nu_j = 3 * (nb - 1) + jcart
|
2004-03-07 21:47:42 +08:00
|
|
|
dy1 (nu_i, nu_j) = dy1 (nu_i, nu_j) + facg * (q (icart) + &
|
|
|
|
g (icart, ng) ) * (q (jcart) + g (jcart, ng) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! the second term
|
|
|
|
!
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
gt2 = gg (ng) * tpiba2
|
2004-03-07 21:47:42 +08:00
|
|
|
if (abs (gt2) > 1.d-8) then
|
|
|
|
fac = - e2 * fpi * tpiba2 / omega * exp ( - gt2 / alpha / 4.d0) / gt2
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
|
|
|
fac = 0.d0
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
nta = ityp (na)
|
|
|
|
fnat = (0.d0, 0.d0)
|
|
|
|
do nb = 1, nat
|
|
|
|
ntb = ityp (nb)
|
2003-01-20 05:58:50 +08:00
|
|
|
arg = tpi * ( (g (1, ng) ) * (tau (1, na) - tau (1, nb) ) + &
|
2004-03-07 21:47:42 +08:00
|
|
|
(g (2, ng) ) * (tau (2, na) - tau (2, nb) ) + &
|
|
|
|
(g (3, ng) ) * (tau (3, na) - tau (3, nb) ) )
|
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
|
|
|
facg = fac * zv (nta) * zv (ntb) * CMPLX (cos (arg), 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
fnat = fnat + facg
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do icart = 1, 3
|
|
|
|
nu_i = 3 * (na - 1) + icart
|
|
|
|
do jcart = 1, 3
|
|
|
|
nu_j = 3 * (na - 1) + jcart
|
2003-01-20 05:58:50 +08:00
|
|
|
dy2 (nu_i, nu_j) = dy2 (nu_i, nu_j) + fnat * g (icart, ng) &
|
|
|
|
* g (jcart, ng)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do nu_i = 1, nmodes
|
|
|
|
do nu_j = 1, nmodes
|
2004-03-07 21:47:42 +08:00
|
|
|
dy3 (nu_i, nu_j) = dy3 (nu_i, nu_j) + dy1 (nu_i, nu_j) - &
|
|
|
|
dy2 (nu_i, nu_j)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! Then there is also a part in real space which is computed here.
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2003-01-20 05:58:50 +08:00
|
|
|
! ... only by the node that contains G=0
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
if (gg (1) > 1.d-8) goto 100
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
rmax = 5.d0 / sqrt (alpha) / alat
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! with this choice terms up to ZiZj*erfc(5) are counted (erfc(5)=2x10^-1
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
nta = ityp (na)
|
|
|
|
do nb = 1, nat
|
|
|
|
ntb = ityp (nb)
|
|
|
|
do icart = 1, 3
|
|
|
|
dtau (icart) = tau (icart, na) - tau (icart, nb)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! generates nearest-neighbors shells r(i)=R(i)-dtau(i)
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call rgen (dtau, rmax, mxr, at, bg, r, r2, nrm)
|
|
|
|
do nr = 1, nrm
|
|
|
|
rr = sqrt (r2 (nr) ) * alat
|
|
|
|
ar = sqrt (alpha) * rr
|
2004-03-07 21:47:42 +08:00
|
|
|
qrg = tpi * (q (1) * (r (1, nr) + dtau (1) ) + &
|
|
|
|
q (2) * (r (2, nr) + dtau (2) ) + &
|
|
|
|
q (3) * (r (3, nr) + dtau (3) ) )
|
|
|
|
d2f = (3.d0 * erfc (ar) + sqrt (8.d0 / tpi) * ar * &
|
|
|
|
(3.d0 + 2.d0 * ar**2) * exp ( - ar**2) ) / rr**5
|
2003-01-20 05:58:50 +08:00
|
|
|
df = ( - erfc (ar) - sqrt (8.d0 / tpi) * ar * exp ( - ar**2) ) &
|
|
|
|
/ rr**3
|
2003-02-08 00:04:36 +08:00
|
|
|
do icart = 1, 3
|
|
|
|
na_icart = 3 * (na - 1) + icart
|
|
|
|
nb_icart = 3 * (nb - 1) + icart
|
|
|
|
do jcart = 1, 3
|
|
|
|
nb_jcart = 3 * (nb - 1) + jcart
|
|
|
|
na_jcart = 3 * (na - 1) + jcart
|
2004-03-07 21:47:42 +08:00
|
|
|
dy3 (na_icart, nb_jcart) = dy3 (na_icart, nb_jcart) + &
|
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
|
|
|
e2 * zv (nta) * zv (ntb) * CMPLX (cos (qrg), sin (qrg))&
|
2004-03-07 21:47:42 +08:00
|
|
|
* (d2f * alat * r (icart, nr) * alat * r (jcart, nr) )
|
|
|
|
dy3 (na_icart, na_jcart) = dy3 (na_icart, na_jcart) - &
|
|
|
|
e2 * zv (nta) * zv (ntb) * (d2f * alat * r (icart, nr) *&
|
|
|
|
alat * r (jcart, nr) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
dy3 (na_icart, nb_icart) = dy3 (na_icart, nb_icart) + e2 * &
|
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
|
|
|
zv (nta) * zv (ntb) * CMPLX (cos (qrg), sin (qrg) ) * df
|
2004-03-07 21:47:42 +08:00
|
|
|
dy3 (na_icart, na_icart) = dy3 (na_icart, na_icart) - e2 * &
|
|
|
|
zv (nta) * zv (ntb) * df
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2003-02-08 00:04:36 +08:00
|
|
|
100 continue
|
2008-04-19 01:05:56 +08:00
|
|
|
call mp_sum ( dy3, intra_pool_comm )
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
|
|
|
!
|
|
|
|
! The dynamical matrix was computed in cartesian axis and now we put
|
|
|
|
! it on the basis of the modes
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do nu_i = 1, nmodes
|
|
|
|
do nu_j = 1, nmodes
|
|
|
|
work = (0.d0, 0.d0)
|
|
|
|
do nb_jcart = 1, 3 * nat
|
|
|
|
do na_icart = 1, 3 * nat
|
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
|
|
|
work = work + CONJG(u (na_icart, nu_i) ) * &
|
2004-03-07 21:47:42 +08:00
|
|
|
dy3 (na_icart, nb_jcart) * u (nb_jcart, nu_j)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
dyn (nu_i, nu_j) = dyn (nu_i, nu_j) - work
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
call stop_clock ('d2ionq')
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine d2ionq
|