quantum-espresso/PH/sym_def.f90

73 lines
1.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 sym_def (def, irr)
!---------------------------------------------------------------------
! Symmetrizes the first order changes of the Fermi energies of an
! irreducible representation. These objects are defined complex because
! perturbations may be complex
!
! Used in the q=0 metallic case only.
!
USE kinds, only : DP
use phcom
implicit none
integer :: irr
! input: the representation under consideration
complex(DP) :: def (npertx)
! inp/out: the fermi energy changes
integer :: ipert, jpert, isym, irot
! counter on perturbations
! counter on perturbations
! counter on symmetries
! the rotation
complex(DP) :: w_def (npertx)
! the fermi energy changes (work array)
if (nsymq == 1 .and. (.not.minus_q) ) return
!
! first the symmetrization S(irotmq)*q = -q + Gi if necessary
!
if (minus_q) then
w_def = (0.d0, 0.d0)
do ipert = 1, npert (irr)
do jpert = 1, npert (irr)
w_def (ipert) = w_def (ipert) + tmq (jpert, ipert, irr) &
* def (jpert)
enddo
enddo
do ipert = 1, npert (irr)
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
def (ipert) = 0.5d0 * (def (ipert) + CONJG(w_def (ipert) ) )
enddo
endif
!
! Here we symmetrize with respect to the small group of q
!
w_def = (0.d0, 0.d0)
do ipert = 1, npert (irr)
do isym = 1, nsymq
irot = irgq (isym)
do jpert = 1, npert (irr)
w_def (ipert) = w_def (ipert) + t (jpert, ipert, irot, irr) &
* def (jpert)
enddo
enddo
enddo
!
! normalize and exit
!
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
def = w_def / DBLE(nsymq)
return
end subroutine sym_def