pw90common_fourier_R_to_k_vec Subroutine

public subroutine pw90common_fourier_R_to_k_vec(kpt, OO_R, OO_true, OO_pseudo)

Uses

  • proc~~pw90common_fourier_r_to_k_vec~~UsesGraph proc~pw90common_fourier_r_to_k_vec pw90common_fourier_R_to_k_vec module~w90_ws_distance w90_ws_distance proc~pw90common_fourier_r_to_k_vec->module~w90_ws_distance module~w90_parameters w90_parameters proc~pw90common_fourier_r_to_k_vec->module~w90_parameters module~w90_constants w90_constants proc~pw90common_fourier_r_to_k_vec->module~w90_constants module~w90_ws_distance->module~w90_parameters module~w90_ws_distance->module~w90_constants module~w90_parameters->module~w90_constants module~w90_io w90_io module~w90_parameters->module~w90_io module~w90_io->module~w90_constants

For OO_true (true vector):

Arguments

Type IntentOptional AttributesName
real(kind=dp) :: kpt(3)
complex(kind=dp), intent(in), dimension(:, :, :, :):: OO_R
complex(kind=dp), intent(out), optional dimension(:, :, :):: OO_true
complex(kind=dp), intent(out), optional dimension(:, :, :):: OO_pseudo

Called by

proc~~pw90common_fourier_r_to_k_vec~~CalledByGraph proc~pw90common_fourier_r_to_k_vec pw90common_fourier_R_to_k_vec proc~berry_get_kubo_k berry_get_kubo_k proc~berry_get_kubo_k->proc~pw90common_fourier_r_to_k_vec proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->proc~pw90common_fourier_r_to_k_vec proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~gyrotropic_get_k_list->proc~berry_get_imfgh_klist proc~berry_get_imf_klist berry_get_imf_klist proc~gyrotropic_get_k_list->proc~berry_get_imf_klist proc~berry_get_shc_klist berry_get_shc_klist proc~berry_get_shc_klist->proc~pw90common_fourier_r_to_k_vec proc~berry_get_imfgh_klist->proc~pw90common_fourier_r_to_k_vec proc~berry_get_imf_klist->proc~berry_get_imfgh_klist proc~k_slice k_slice proc~k_slice->proc~berry_get_shc_klist proc~k_slice->proc~berry_get_imfgh_klist proc~k_slice->proc~berry_get_imf_klist proc~k_path k_path proc~k_path->proc~berry_get_shc_klist proc~k_path->proc~berry_get_imfgh_klist proc~k_path->proc~berry_get_imf_klist proc~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~gyrotropic_get_k_list proc~berry_main berry_main proc~berry_main->proc~berry_get_kubo_k proc~berry_main->proc~berry_get_shc_klist proc~berry_main->proc~berry_get_imfgh_klist proc~berry_main->proc~berry_get_imf_klist

Contents


Source Code

  subroutine pw90common_fourier_R_to_k_vec(kpt, OO_R, OO_true, OO_pseudo)
    !====================================================================!
    !                                                                    !
    !! For OO_true (true vector):
    !! $${\vec O}_{ij}(k) = \sum_R e^{+ik.R} {\vec O}_{ij}(R)$$
    !                                                                    !
    !====================================================================!

    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), optional, dimension(:, :, :), intent(out)   :: OO_true
    complex(kind=dp), optional, dimension(:, :, :), intent(out)   :: OO_pseudo

    integer          :: ir, i, j, ideg
    real(kind=dp)    :: rdotk
    complex(kind=dp) :: phase_fac

    if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
    if (present(OO_true)) OO_true = cmplx_0
    if (present(OO_pseudo)) OO_pseudo = 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 (present(OO_true)) then
              OO_true(i, j, 1) = OO_true(i, j, 1) + phase_fac*OO_R(i, j, ir, 1)
              OO_true(i, j, 2) = OO_true(i, j, 2) + phase_fac*OO_R(i, j, ir, 2)
              OO_true(i, j, 3) = OO_true(i, j, 3) + phase_fac*OO_R(i, j, ir, 3)
            endif
            if (present(OO_pseudo)) then
              OO_pseudo(i, j, 1) = OO_pseudo(i, j, 1) &
                                   + cmplx_i*crdist_ws(2, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 3) &
                                   - cmplx_i*crdist_ws(3, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 2)
              OO_pseudo(i, j, 2) = OO_pseudo(i, j, 2) &
                                   + cmplx_i*crdist_ws(3, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 1) &
                                   - cmplx_i*crdist_ws(1, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 3)
              OO_pseudo(i, j, 3) = OO_pseudo(i, j, 3) &
                                   + cmplx_i*crdist_ws(1, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 2) &
                                   - cmplx_i*crdist_ws(2, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir, 1)
            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 (present(OO_true)) then
          OO_true(:, :, 1) = OO_true(:, :, 1) + phase_fac*OO_R(:, :, ir, 1)
          OO_true(:, :, 2) = OO_true(:, :, 2) + phase_fac*OO_R(:, :, ir, 2)
          OO_true(:, :, 3) = OO_true(:, :, 3) + phase_fac*OO_R(:, :, ir, 3)
        endif
        if (present(OO_pseudo)) then
          OO_pseudo(:, :, 1) = OO_pseudo(:, :, 1) &
                               + cmplx_i*crvec(2, ir)*phase_fac*OO_R(:, :, ir, 3) &
                               - cmplx_i*crvec(3, ir)*phase_fac*OO_R(:, :, ir, 2)
          OO_pseudo(:, :, 2) = OO_pseudo(:, :, 2) &
                               + cmplx_i*crvec(3, ir)*phase_fac*OO_R(:, :, ir, 1) &
                               - cmplx_i*crvec(1, ir)*phase_fac*OO_R(:, :, ir, 3)
          OO_pseudo(:, :, 3) = OO_pseudo(:, :, 3) &
                               + cmplx_i*crvec(1, ir)*phase_fac*OO_R(:, :, ir, 2) &
                               - cmplx_i*crvec(2, ir)*phase_fac*OO_R(:, :, ir, 1)
        endif
      endif
    enddo

  end subroutine pw90common_fourier_R_to_k_vec