quantum-espresso/LR_Modules/mode_group.f90

118 lines
3.7 KiB
Fortran
Raw Normal View History

!
! Copyright (C) 2001-2008 Quantum ESPRESSO 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 mode_group &
(modenum, xq, at, bg, nat, nrot, s, irt, minus_q, rtau, sym)
!-----------------------------------------------------------------------
!
! This routine selects, among the symmetry matrices of the point group
! of a crystal, the symmetry operations which leave a given mode unchanged
! For the moment it assumes that the mode modenum displaces the atom
! modenum/3 in the direction mod(modenum,3)+1
!
USE kinds, ONLY : DP
USE constants, ONLY : tpi
implicit none
integer, intent(in) :: nat, s (3, 3, 48), irt (48, nat), nrot, modenum
! nat : the number of atoms of the system
! s : the symmetry matrices
! irt : the rotated atom
! nrot: number of symmetry operations
! modenum: which displacement pattern
real(DP), intent(in) :: xq (3), rtau (3, 48, nat), bg (3, 3), at (3, 3)
! xq : the q point
! rtau: the translations of each atom
! bg : the reciprocal lattice vectors
! at : the direct lattice vectors
logical, intent(in) :: minus_q
! if true Sq=>-q+G symmetry is used
logical, intent(inout) :: sym (48)
! on input: .true. if symm. op. has to be tested
! on output: .true. if symm. op. does not change mode modenum
!
integer :: isym, nas, ipols, na, sna, ipol, jpol
! counters
real(DP) :: arg
! auxiliary
complex(DP), allocatable :: u (:,:)
! the original pattern
complex(DP) :: fase, sum
! the phase of the mode
! check for orthogonality
complex(DP), allocatable :: work_u (:,:), work_ru (:,:)
! the working pattern
! the rotated working pattern
allocate(u(3, nat), work_u(3, nat), work_ru (3, nat))
if (modenum > 3*nat .or. modenum < 1) call errore ('mode_group', &
'wrong modenum', 1)
nas = (modenum - 1) / 3 + 1
ipols = mod (modenum - 1, 3) + 1
u (:,:) = (0.d0, 0.d0)
u (ipols, nas) = (1.d0, 0.d0)
do na = 1, nat
call trnvecc (u (1, na), at, bg, - 1)
enddo
do isym = 1, nrot
if (sym (isym) ) then
do na = 1, nat
do ipol = 1, 3
work_u (ipol, na) = u (ipol, na)
enddo
enddo
work_ru (:,:) = (0.d0, 0.d0)
do na = 1, nat
sna = irt (isym, na)
arg = 0.d0
do ipol = 1, 3
arg = arg + xq (ipol) * rtau (ipol, isym, na)
enddo
arg = arg * tpi
if (isym.eq.nrot.and.minus_q) then
fase = CMPLX(cos (arg), sin (arg) ,kind=DP)
else
fase = CMPLX(cos (arg), - sin (arg) ,kind=DP)
endif
do ipol = 1, 3
do jpol = 1, 3
work_ru (ipol, sna) = work_ru (ipol, sna) + s (jpol, ipol, &
isym) * work_u (jpol, na) * fase
enddo
enddo
enddo
!
! Transform back the rotated pattern
!
do na = 1, nat
call trnvecc (work_ru (1, na), at, bg, 1)
call trnvecc (work_u (1, na), at, bg, 1)
enddo
!
! only if the pattern remain the same up to a phase we keep
! the symmetry
!
sum = (0.d0, 0.d0)
do na = 1, nat
do ipol = 1, 3
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
sum = sum + CONJG(work_u (ipol, na) ) * work_ru (ipol, na)
enddo
enddo
sum = abs (sum)
if (abs (sum - 1.d0) .gt.1.d-7) sym (isym) = .false.
endif
enddo
deallocate ( work_ru, work_u, u)
return
end subroutine mode_group