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 symdynph_gq (xq, phi, s, invs, rtau, irt, irgq, nsymq, &
|
|
|
|
nat, irotmq, minus_q)
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine receives as input an unsymmetrized dynamical
|
|
|
|
! matrix expressed on the crystal axes and imposes the symmetry
|
|
|
|
! of the small group of q. Furthermore it imposes also the symmetry
|
|
|
|
! q -> -q+G if present.
|
|
|
|
!
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2004-03-07 21:47:42 +08:00
|
|
|
USE constants, ONLY: tpi
|
2003-02-08 00:04:36 +08:00
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! The dummy variables
|
|
|
|
!
|
|
|
|
integer :: nat, s (3, 3, 48), irt (48, nat), irgq (48), invs (48), &
|
|
|
|
nsymq, irotmq
|
|
|
|
! input: the number of atoms
|
|
|
|
! input: the symmetry matrices
|
|
|
|
! input: the rotated of each vector
|
|
|
|
! input: the small group of q
|
|
|
|
! input: the inverse of each matrix
|
|
|
|
! input: the order of the small gro
|
|
|
|
! input: the rotation sending q ->
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: xq (3), rtau (3, 48, nat)
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the q point
|
|
|
|
! input: the R associated at each t
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
logical :: minus_q
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: true if a symmetry q->-q+G
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: phi (3, 3, nat, nat)
|
2003-01-20 05:58:50 +08:00
|
|
|
! inp/out: the matrix to symmetrize
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
! local variables
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
integer :: isymq, sna, snb, irot, na, nb, ipol, jpol, lpol, kpol, &
|
|
|
|
iflb (nat, nat)
|
2004-03-07 21:47:42 +08:00
|
|
|
! counters, indices, work space
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: arg
|
2003-01-20 05:58:50 +08:00
|
|
|
! the argument of the phase
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: phip (3, 3, nat, nat), work (3, 3), fase, faseq (48)
|
2004-03-07 21:47:42 +08:00
|
|
|
! work space, phase factors
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! We start by imposing hermiticity
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
do nb = 1, nat
|
|
|
|
do ipol = 1, 3
|
|
|
|
do jpol = 1, 3
|
2003-01-20 05:58:50 +08:00
|
|
|
phi (ipol, jpol, na, nb) = 0.5d0 * (phi (ipol, jpol, na, 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
|
|
|
+ CONJG(phi (jpol, ipol, nb, na) ) )
|
|
|
|
phi (jpol, ipol, nb, na) = CONJG(phi (ipol, jpol, na, nb) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! If no other symmetry is present we quit here
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
if ( (nsymq == 1) .and. (.not.minus_q) ) return
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Then we impose the symmetry q -> -q+G if present
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
if (minus_q) then
|
|
|
|
do na = 1, nat
|
|
|
|
do nb = 1, nat
|
|
|
|
do ipol = 1, 3
|
|
|
|
do jpol = 1, 3
|
2004-03-07 21:47:42 +08:00
|
|
|
work(:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
sna = irt (irotmq, na)
|
|
|
|
snb = irt (irotmq, nb)
|
|
|
|
arg = 0.d0
|
|
|
|
do kpol = 1, 3
|
2004-03-07 21:47:42 +08:00
|
|
|
arg = arg + (xq (kpol) * (rtau (kpol, irotmq, na) - &
|
|
|
|
rtau (kpol, irotmq, nb) ) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
arg = arg * tpi
|
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) )
|
2005-01-07 20:57:41 +08:00
|
|
|
#if defined __ALTIX
|
|
|
|
!DIR$ unroll (0)
|
|
|
|
#endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do kpol = 1, 3
|
|
|
|
do lpol = 1, 3
|
2004-03-07 21:47:42 +08:00
|
|
|
work (ipol, jpol) = work (ipol, jpol) + &
|
|
|
|
s (ipol, kpol, irotmq) * s (jpol, lpol, irotmq) &
|
|
|
|
* phi (kpol, lpol, sna, snb) * fase
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
phip (ipol, jpol, na, nb) = (phi (ipol, jpol, na, 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
|
|
|
CONJG( work (ipol, jpol) ) ) * 0.5d0
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
phi = phip
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
|
|
|
|
!
|
|
|
|
! Here we symmetrize with respect to the small group of q
|
|
|
|
!
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nsymq == 1) return
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
iflb (:, :) = 0
|
2003-02-08 00:04:36 +08:00
|
|
|
do na = 1, nat
|
|
|
|
do nb = 1, nat
|
2004-03-07 21:47:42 +08:00
|
|
|
if (iflb (na, nb) == 0) then
|
|
|
|
work(:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
do isymq = 1, nsymq
|
|
|
|
irot = irgq (isymq)
|
|
|
|
sna = irt (irot, na)
|
|
|
|
snb = irt (irot, nb)
|
|
|
|
arg = 0.d0
|
|
|
|
do ipol = 1, 3
|
2004-03-07 21:47:42 +08:00
|
|
|
arg = arg + (xq (ipol) * (rtau (ipol, irot, na) - &
|
|
|
|
rtau (ipol, irot, nb) ) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
arg = arg * tpi
|
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
|
|
|
faseq (isymq) = CMPLX (cos (arg), sin (arg) )
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipol = 1, 3
|
|
|
|
do jpol = 1, 3
|
|
|
|
do kpol = 1, 3
|
|
|
|
do lpol = 1, 3
|
2004-03-07 21:47:42 +08:00
|
|
|
work (ipol, jpol) = work (ipol, jpol) + &
|
|
|
|
s (ipol, kpol, irot) * s (jpol, lpol, irot) &
|
|
|
|
* phi (kpol, lpol, sna, snb) * faseq (isymq)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do isymq = 1, nsymq
|
|
|
|
irot = irgq (isymq)
|
|
|
|
sna = irt (irot, na)
|
|
|
|
snb = irt (irot, nb)
|
|
|
|
do ipol = 1, 3
|
|
|
|
do jpol = 1, 3
|
|
|
|
phi (ipol, jpol, sna, snb) = (0.d0, 0.d0)
|
|
|
|
do kpol = 1, 3
|
|
|
|
do lpol = 1, 3
|
2003-01-20 05:58:50 +08:00
|
|
|
phi (ipol, jpol, sna, snb) = phi (ipol, jpol, sna, snb) &
|
|
|
|
+ s (ipol, kpol, invs (irot) ) * s (jpol, lpol, invs (irot) ) &
|
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 (kpol, lpol) * CONJG(faseq (isymq) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
iflb (sna, snb) = 1
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
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
|
|
|
phi (:, :, :, :) = phi (:, :, :, :) / DBLE(nsymq)
|
2003-02-08 00:04:36 +08:00
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine symdynph_gq
|