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 .
|
|
|
|
!
|
2004-10-26 17:32:48 +08:00
|
|
|
#include "f_defs.h"
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
|
|
|
|
subroutine mode_group (modenum, xq, at, bg, nat, nrot, s, irt, &
|
|
|
|
rtau, sym, minus_q)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine selects, among the symmetry matrices of the point group
|
|
|
|
! of a crystal, the symmetry operations which leave a given mode unchang
|
|
|
|
! For the moment it assume that the mode modenum displaces the atom
|
|
|
|
! modenum/3 in the direction mod(modenum,3)+1
|
|
|
|
! Also the minus_q operation is tested.
|
|
|
|
!
|
|
|
|
! input-output variables
|
|
|
|
!
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds
|
2003-02-08 00:04:36 +08:00
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: nat, s (3, 3, 48), irt (48, nat), nrot, modenum
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the number of atoms of the system
|
|
|
|
! input: the symmetry matrices
|
|
|
|
! input: the rotated atom
|
|
|
|
! input: number of symmetry operations
|
|
|
|
! input: the displacement pattern
|
|
|
|
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: xq (3), rtau (3, 48, nat), bg (3, 3), at (3, 3)
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the q point
|
|
|
|
! input: the translations of each atom
|
|
|
|
! input: the reciprocal lattice vectors
|
|
|
|
! input: the direct lattice vectors
|
2003-02-08 00:04:36 +08:00
|
|
|
logical :: minus_q, sym (48)
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: if true minus_q symmetry is used
|
|
|
|
! input-output: .true. if symm. op. do not change
|
|
|
|
! mode
|
|
|
|
!
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: isym, nas, ipols, na, sna, ipol, jpol
|
2003-01-20 05:58:50 +08:00
|
|
|
! counters
|
|
|
|
! counter on polarizations
|
|
|
|
! counter on polarizations
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP), parameter :: tpi = 2.0d0 * 3.14159265358979d0
|
|
|
|
real(DP) :: arg
|
2003-01-20 05:58:50 +08:00
|
|
|
! auxiliary
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: u (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the original pattern
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: fase, sum
|
2003-01-20 05:58:50 +08:00
|
|
|
! the phase of the mode
|
|
|
|
! check for orthogonality
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: work_u (:,:), work_ru (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the working pattern
|
|
|
|
! the rotated working pattern
|
|
|
|
|
|
|
|
|
|
|
|
allocate(u(3, nat), work_u(3, nat), work_ru (3, nat))
|
|
|
|
|
2003-02-21 22:57:00 +08:00
|
|
|
if (modenum.gt.3 * nat.or.modenum.lt.1) call errore ('mode_group', &
|
2003-01-20 05:58:50 +08:00
|
|
|
'wrong modenum', 1)
|
2003-02-08 00:04:36 +08:00
|
|
|
nas = (modenum - 1) / 3 + 1
|
|
|
|
ipols = mod (modenum - 1, 3) + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
u (:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
u (ipols, nas) = (1.d0, 0.d0)
|
|
|
|
do na = 1, nat
|
|
|
|
call trnvecc (u (1, na), at, bg, - 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nrot
|
|
|
|
if (sym (isym) ) then
|
|
|
|
do na = 1, nat
|
|
|
|
do ipol = 1, 3
|
|
|
|
work_u (ipol, na) = u (ipol, na)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
work_ru (:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
sna = irt (isym, na)
|
|
|
|
arg = 0.d0
|
|
|
|
do ipol = 1, 3
|
|
|
|
arg = arg + xq (ipol) * rtau (ipol, isym, na)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
arg = arg * tpi
|
|
|
|
if (isym.eq.nrot.and.minus_q) then
|
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
|
|
|
fase = CMPLX (cos (arg), sin (arg) )
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
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
|
|
|
fase = CMPLX (cos (arg), - sin (arg) )
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipol = 1, 3
|
|
|
|
do jpol = 1, 3
|
2003-01-20 05:58:50 +08:00
|
|
|
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
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
call trnvecc (work_ru (1, na), at, bg, 1)
|
|
|
|
call trnvecc (work_u (1, na), at, bg, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! only if the pattern remain the same ap to a phase we keep
|
|
|
|
! the symmetry
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
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)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
sum = abs (sum)
|
|
|
|
if (abs (sum - 1.d0) .gt.1.d-7) sym (isym) = .false.
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate ( work_ru, work_u, u)
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
end subroutine mode_group
|
|
|
|
|