For alpha=0: O_ij(R) --> O_ij(k) = sum_R e^{+ik.R}*O_ij(R)
For alpha=1,2,3: sum_R [cmplx_iR_alphae^{+ik.R}*O_ij(R)] where R_alpha is a Cartesian component of R REMOVE EVENTUALLY (replace with pw90common_fourier_R_to_k_new)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp) | :: | kpt(3) | ||||
complex(kind=dp), | intent(in), | dimension(:, :, :) | :: | OO_R | ||
complex(kind=dp), | intent(out), | dimension(:, :) | :: | OO | ||
integer | :: | alpha |
subroutine pw90common_fourier_R_to_k(kpt, OO_R, OO, alpha)
!=========================================================!
! !
!! For alpha=0:
!! O_ij(R) --> O_ij(k) = sum_R e^{+ik.R}*O_ij(R)
!!
!! For alpha=1,2,3:
!! sum_R [cmplx_i*R_alpha*e^{+ik.R}*O_ij(R)]
!! where R_alpha is a Cartesian component of R
!! ***REMOVE EVENTUALLY*** (replace with pw90common_fourier_R_to_k_new)
! !
!=========================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_parameters, only: num_kpts, kpt_latt, num_wann, use_ws_distance
use w90_ws_distance, only: irdist_ws, crdist_ws, &
wdist_ndeg, ws_translate_dist
implicit none
! Arguments
!
real(kind=dp) :: kpt(3)
complex(kind=dp), dimension(:, :, :), intent(in) :: OO_R
complex(kind=dp), dimension(:, :), intent(out) :: OO
integer :: alpha
integer :: ir, i, j, ideg
real(kind=dp) :: rdotk
complex(kind=dp) :: phase_fac
if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
OO(:, :) = cmplx_0
do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, wdist_ndeg(i, j, ir)
rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir), dp))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
if (alpha == 0) then
OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
elseif (alpha == 1 .or. alpha == 2 .or. alpha == 3) then
OO(i, j) = OO(i, j) + cmplx_i*crdist_ws(alpha, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir)
else
stop 'wrong value of alpha in pw90common_fourier_R_to_k'
endif
enddo
enddo
enddo
else
! [lp] Original code, without IJ-dependent shift:
rdotk = twopi*dot_product(kpt(:), irvec(:, ir))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir), dp)
if (alpha == 0) then
OO(:, :) = OO(:, :) + phase_fac*OO_R(:, :, ir)
elseif (alpha == 1 .or. alpha == 2 .or. alpha == 3) then
OO(:, :) = OO(:, :) + &
cmplx_i*crvec(alpha, ir)*phase_fac*OO_R(:, :, ir)
else
stop 'wrong value of alpha in pw90common_fourier_R_to_k'
endif
endif
enddo
end subroutine pw90common_fourier_R_to_k