quantum-espresso/D3/d3_exc.f90

98 lines
2.9 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 d3_exc
!-----------------------------------------------------------------------
!
! Calculates the contribution to the derivative of the dynamical
! matrix due to the third derivative of the exchange and correlation
! energy
!
USE ions_base, ONLY : nat
USE kinds, ONLY : DP
USE pwcom
USE scf, only : rho, rho_core
USE phcom
USE d3com
USE io_global, ONLY : ionode_id
USE mp_global, ONLY : inter_pool_comm, my_image_id, me_pool, &
root_image, npool, intra_pool_comm
USE mp, ONLY : mp_bcast, mp_sum
IMPLICIT NONE
INTEGER :: errcode, ir, ipert, jpert, kpert, npert1, npert2
REAL (DP) :: d2mxc, rhotot, xq0 (3)
REAL (DP), ALLOCATABLE :: d2muxc (:)
COMPLEX (DP) :: aux
COMPLEX (DP), ALLOCATABLE :: work1 (:), work2 (:), &
work3 (:), d3dyn1 (:,:,:)
ALLOCATE (d2muxc( nrxx))
ALLOCATE (work1 ( nrxx))
ALLOCATE (work2 ( nrxx))
ALLOCATE (work3 ( nrxx))
ALLOCATE (d3dyn1( 3*nat, 3*nat, 3*nat))
IF ( me_pool == root_image ) THEN
!
! Calculates third derivative of Exc
!
d2muxc(:) = 0.d0
DO ir = 1, nrxx
rhotot = rho%of_r (ir, 1) + rho_core (ir)
IF (rhotot > 1.d-30) d2muxc (ir) = d2mxc (rhotot)
IF (rhotot < - 1.d-30) d2muxc (ir) = - d2mxc ( - rhotot)
ENDDO
!
! Calculates the contribution to d3dyn
!
d3dyn1 (:,:,:) = (0.d0, 0.d0)
DO ipert = 1, 3 * nat
IF (q0mode (ipert) ) THEN
CALL davcio_drho (work1, lrdrho, iud0rho, ipert, - 1)
DO jpert = 1, 3 * nat
CALL davcio_drho (work2, lrdrho, iudrho, jpert, - 1)
DO kpert = 1, 3 * nat
CALL davcio_drho (work3, lrdrho, iudrho, kpert, - 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
aux = CMPLX (0.d0, 0.d0)
DO ir = 1, nrxx
aux = aux + &
d2muxc (ir) * work1 (ir) * &
CONJG (work2 (ir) ) * work3 (ir)
ENDDO
!
CALL mp_sum ( aux, intra_pool_comm )
!
d3dyn1 (ipert, jpert, kpert) = omega * aux / (nr1 * nr2 * nr3)
!
ENDDO
ENDDO
ENDIF
ENDDO
!
END IF
!
IF ( npool /= 1 ) CALL mp_bcast( d3dyn1, ionode_id, inter_pool_comm )
!
d3dyn = d3dyn + d3dyn1
d3dyn_aux9 = d3dyn1
!
DEALLOCATE (d2muxc)
DEALLOCATE (work1)
DEALLOCATE (work2)
DEALLOCATE (work3)
DEALLOCATE (d3dyn1)
!
RETURN
!
END SUBROUTINE d3_exc