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
|
|
|
!
|
|
|
|
!---------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine symdvscf (nper, irr, dvtosym)
|
2003-01-20 05:58:50 +08:00
|
|
|
!---------------------------------------------------------------------
|
|
|
|
! symmetrize the self-consistent potential of the perturbations
|
|
|
|
! belonging to an irreproducible representation
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only : DP
|
2003-02-08 00:04:36 +08:00
|
|
|
use phcom
|
|
|
|
implicit none
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: nper, irr
|
2003-01-20 05:58:50 +08:00
|
|
|
! the number of perturbations
|
|
|
|
! the representation under conside
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: dvtosym (nrx1, nrx2, nrx3, nspin, nper)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the potential to symmetriz
|
|
|
|
|
|
|
|
integer :: is, ri, rj, rk, i, j, k, ipert, jpert, ipol, isym, &
|
|
|
|
irot
|
|
|
|
! counter on spin polarizations
|
|
|
|
!
|
|
|
|
! the rotated points
|
|
|
|
!
|
|
|
|
!
|
|
|
|
! counter on mesh points
|
|
|
|
!
|
|
|
|
! counter on perturbations
|
|
|
|
! counter on perturbations
|
|
|
|
! counter on polarizations
|
|
|
|
! counter on symmetries
|
|
|
|
! the rotation
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: g1 (48), g2 (48), g3 (48), in1, in2, in3
|
2003-01-20 05:58:50 +08:00
|
|
|
! used to construct the phases
|
|
|
|
! auxiliary variables
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP), allocatable :: dvsym (:,:,:,:)
|
2004-03-07 21:47:42 +08:00
|
|
|
! the symmetrized potential
|
2005-08-28 22:09:42 +08:00
|
|
|
complex(DP) :: aux2, term (3, 48), phase (48)
|
2003-01-20 05:58:50 +08:00
|
|
|
! auxiliary space
|
|
|
|
! the multiplication factor
|
|
|
|
! the phase factor
|
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (nsymq == 1.and. (.not.minus_q) ) return
|
2003-02-08 00:04:36 +08:00
|
|
|
call start_clock ('symdvscf')
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (dvsym( nrx1 , nrx2 , nrx3 , nper))
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! if necessary we symmetrize with respect to S(irotmq)*q = -q + Gi
|
|
|
|
!
|
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
|
|
|
in1 = tpi / DBLE (nr1)
|
|
|
|
in2 = tpi / DBLE (nr2)
|
|
|
|
in3 = tpi / DBLE (nr3)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (minus_q) then
|
|
|
|
g1 (1) = 0.d0
|
|
|
|
g2 (1) = 0.d0
|
|
|
|
g3 (1) = 0.d0
|
|
|
|
do ipol = 1, 3
|
|
|
|
g1 (1) = g1 (1) + gimq (ipol) * in1 * at (ipol, 1)
|
|
|
|
g2 (1) = g2 (1) + gimq (ipol) * in2 * at (ipol, 2)
|
|
|
|
g3 (1) = g3 (1) + gimq (ipol) * in3 * at (ipol, 3)
|
2003-01-20 05:58:50 +08:00
|
|
|
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
|
|
|
term (1, 1) = CMPLX (cos (g1 (1) ), sin (g1 (1) ) )
|
|
|
|
term (2, 1) = CMPLX (cos (g2 (1) ), sin (g2 (1) ) )
|
|
|
|
term (3, 1) = CMPLX (cos (g3 (1) ), sin (g3 (1) ) )
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
|
|
|
phase (1) = (1.d0, 0.d0)
|
|
|
|
do k = 1, nr3
|
|
|
|
do j = 1, nr2
|
|
|
|
do i = 1, nr1
|
2003-01-20 05:58:50 +08:00
|
|
|
ri = s (1, 1, irotmq) * (i - 1) + s (2, 1, irotmq) * (j - 1) &
|
|
|
|
+ s (3, 1, irotmq) * (k - 1) - ftau (1, irotmq)
|
2003-02-08 00:04:36 +08:00
|
|
|
ri = mod (ri, nr1) + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
if (ri < 1) ri = ri + nr1
|
2003-01-20 05:58:50 +08:00
|
|
|
rj = s (1, 2, irotmq) * (i - 1) + s (2, 2, irotmq) * (j - 1) &
|
|
|
|
+ s (3, 2, irotmq) * (k - 1) - ftau (2, irotmq)
|
2003-02-08 00:04:36 +08:00
|
|
|
rj = mod (rj, nr2) + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
if (rj < 1) rj = rj + nr2
|
2003-01-20 05:58:50 +08:00
|
|
|
rk = s (1, 3, irotmq) * (i - 1) + s (2, 3, irotmq) * (j - 1) &
|
|
|
|
+ s (3, 3, irotmq) * (k - 1) - ftau (3, irotmq)
|
2003-02-08 00:04:36 +08:00
|
|
|
rk = mod (rk, nr3) + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (rk < 1) rk = rk + nr3
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipert = 1, nper
|
|
|
|
aux2 = (0.d0, 0.d0)
|
|
|
|
do jpert = 1, nper
|
2004-03-07 21:47:42 +08:00
|
|
|
aux2 = aux2 + tmq (jpert, ipert, irr) * &
|
|
|
|
dvtosym (ri, rj, rk, is, jpert) * phase (1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
dvsym (i, j, k, ipert) = (dvtosym (i, j, k, is, ipert) +&
|
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(aux2) ) * 0.5d0
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
phase (1) = phase (1) * term (1, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
phase (1) = phase (1) * term (2, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
phase (1) = phase (1) * term (3, 1)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipert = 1, nper
|
2004-03-07 21:47:42 +08:00
|
|
|
dvtosym(:, :, :, is, ipert) = dvsym (:, :, :, ipert)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
! Here we symmetrize with respect to the small group of q
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nsymq
|
|
|
|
g1 (isym) = 0.d0
|
|
|
|
g2 (isym) = 0.d0
|
|
|
|
g3 (isym) = 0.d0
|
|
|
|
do ipol = 1, 3
|
|
|
|
g1 (isym) = g1 (isym) + gi (ipol, isym) * in1 * at (ipol, 1)
|
|
|
|
g2 (isym) = g2 (isym) + gi (ipol, isym) * in2 * at (ipol, 2)
|
|
|
|
g3 (isym) = g3 (isym) + gi (ipol, isym) * in3 * at (ipol, 3)
|
2003-01-20 05:58:50 +08:00
|
|
|
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
|
|
|
term (1, isym) = CMPLX (cos (g1 (isym) ), sin (g1 (isym) ) )
|
|
|
|
term (2, isym) = CMPLX (cos (g2 (isym) ), sin (g2 (isym) ) )
|
|
|
|
term (3, isym) = CMPLX (cos (g3 (isym) ), sin (g3 (isym) ) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2004-03-07 21:47:42 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
2004-03-07 21:47:42 +08:00
|
|
|
dvsym(:,:,:,:) = (0.d0, 0.d0)
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nsymq
|
|
|
|
phase (isym) = (1.d0, 0.d0)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do k = 1, nr3
|
|
|
|
do j = 1, nr2
|
|
|
|
do i = 1, nr1
|
|
|
|
do isym = 1, nsymq
|
|
|
|
irot = irgq (isym)
|
2004-03-07 21:47:42 +08:00
|
|
|
ri = s (1, 1, irot) * (i - 1) + s (2, 1, irot) * (j - 1) &
|
|
|
|
+ s (3, 1, irot) * (k - 1) - ftau (1, irot)
|
2003-02-08 00:04:36 +08:00
|
|
|
ri = mod (ri, nr1) + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
if (ri < 1) ri = ri + nr1
|
|
|
|
rj = s (1, 2, irot) * (i - 1) + s (2, 2, irot) * (j - 1) &
|
|
|
|
+ s (3, 2, irot) * (k - 1) - ftau (2, irot)
|
2003-02-08 00:04:36 +08:00
|
|
|
rj = mod (rj, nr2) + 1
|
2004-03-07 21:47:42 +08:00
|
|
|
if (rj < 1) rj = rj + nr2
|
|
|
|
rk = s (1, 3, irot) * (i - 1) + s (2, 3, irot) * (j - 1) &
|
|
|
|
+ s (3, 3, irot) * (k - 1) - ftau (3, irot)
|
2003-02-08 00:04:36 +08:00
|
|
|
rk = mod (rk, nr3) + 1
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2004-03-07 21:47:42 +08:00
|
|
|
if (rk < 1) rk = rk + nr3
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipert = 1, nper
|
|
|
|
do jpert = 1, nper
|
2004-03-07 21:47:42 +08:00
|
|
|
dvsym (i, j, k, ipert) = dvsym (i, j, k, ipert) + &
|
|
|
|
t (jpert, ipert, irot, irr) * &
|
|
|
|
dvtosym (ri, rj, rk, is, jpert) * phase (isym)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nsymq
|
|
|
|
phase (isym) = phase (isym) * term (1, isym)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nsymq
|
|
|
|
phase (isym) = phase (isym) * term (2, isym)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do isym = 1, nsymq
|
|
|
|
phase (isym) = phase (isym) * term (3, isym)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
do ipert = 1, nper
|
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
|
|
|
dvtosym(:,:,:,is,ipert) = dvsym(:,:,:,ipert) / DBLE (nsymq)
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (dvsym)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call stop_clock ('symdvscf')
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine symdvscf
|