pw90common_fourier_R_to_k_new Subroutine

public subroutine pw90common_fourier_R_to_k_new(kpt, OO_R, OO, OO_dx, OO_dy, OO_dz)

Uses

  • proc~~pw90common_fourier_r_to_k_new~~UsesGraph proc~pw90common_fourier_r_to_k_new pw90common_fourier_R_to_k_new module~w90_ws_distance w90_ws_distance proc~pw90common_fourier_r_to_k_new->module~w90_ws_distance module~w90_parameters w90_parameters proc~pw90common_fourier_r_to_k_new->module~w90_parameters module~w90_constants w90_constants proc~pw90common_fourier_r_to_k_new->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: For : where R_{x,y,z} are the Cartesian components of R

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

Contents


Source Code

  subroutine pw90common_fourier_R_to_k_new(kpt, OO_R, OO, OO_dx, OO_dy, OO_dz)
    !=======================================================!
    !                                                       !
    !! For OO:
    !! $$O_{ij}(k) = \sum_R e^{+ik.R}.O_{ij}(R)$$
    !! For $$OO_{dx,dy,dz}$$:
    !! $$\sum_R [i.R_{dx,dy,dz}.e^{+ik.R}.O_{ij}(R)]$$
    !! where R_{x,y,z} are the Cartesian components of R
    !                                                       !
    !=======================================================!

    use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
    use w90_parameters, only: timing_level, 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
    complex(kind=dp), optional, dimension(:, :), intent(out)   :: OO_dx
    complex(kind=dp), optional, dimension(:, :), intent(out)   :: OO_dy
    complex(kind=dp), optional, dimension(:, :), intent(out)   :: OO_dz

    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)) OO = cmplx_0
    if (present(OO_dx)) OO_dx = cmplx_0
    if (present(OO_dy)) OO_dy = cmplx_0
    if (present(OO_dz)) OO_dz = 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)) OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
            if (present(OO_dx)) OO_dx(i, j) = OO_dx(i, j) + &
                                              cmplx_i*crdist_ws(1, ideg, i, j, ir)* &
                                              phase_fac*OO_R(i, j, ir)
            if (present(OO_dy)) OO_dy(i, j) = OO_dy(i, j) + &
                                              cmplx_i*crdist_ws(2, ideg, i, j, ir)* &
                                              phase_fac*OO_R(i, j, ir)
            if (present(OO_dz)) OO_dz(i, j) = OO_dz(i, j) + &
                                              cmplx_i*crdist_ws(3, ideg, i, j, ir)* &
                                              phase_fac*OO_R(i, j, ir)
          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)) OO(:, :) = OO(:, :) + phase_fac*OO_R(:, :, ir)
        if (present(OO_dx)) OO_dx(:, :) = OO_dx(:, :) + &
                                          cmplx_i*crvec(1, ir)*phase_fac*OO_R(:, :, ir)
        if (present(OO_dy)) OO_dy(:, :) = OO_dy(:, :) + &
                                          cmplx_i*crvec(2, ir)*phase_fac*OO_R(:, :, ir)
        if (present(OO_dz)) OO_dz(:, :) = OO_dz(:, :) + &
                                          cmplx_i*crvec(3, ir)*phase_fac*OO_R(:, :, ir)
      endif
    enddo

  end subroutine pw90common_fourier_R_to_k_new