fourier_q_to_R Subroutine

private subroutine fourier_q_to_R(op_q, op_R)

Uses

  • proc~~fourier_q_to_r~~UsesGraph proc~fourier_q_to_r fourier_q_to_R module~w90_postw90_common w90_postw90_common proc~fourier_q_to_r->module~w90_postw90_common module~w90_parameters w90_parameters proc~fourier_q_to_r->module~w90_parameters module~w90_constants w90_constants proc~fourier_q_to_r->module~w90_constants module~w90_postw90_common->module~w90_constants module~w90_comms w90_comms module~w90_postw90_common->module~w90_comms module~w90_parameters->module~w90_constants module~w90_io w90_io module~w90_parameters->module~w90_io module~w90_comms->module~w90_constants module~w90_comms->module~w90_io module~w90_io->module~w90_constants

Fourier transforms Wannier-gauge representation of a given operator O from q-space to R-space:

O_ij(q) --> O_ij(R) = (1/N_kpts) sum_q e^{-iqR} O_ij(q)

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :, :):: op_q

Operator in q-space

complex(kind=dp), intent(out), dimension(:, :, :):: op_R

Operator in R-space


Called by

proc~~fourier_q_to_r~~CalledByGraph proc~fourier_q_to_r fourier_q_to_R proc~get_hh_r get_HH_R proc~get_hh_r->proc~fourier_q_to_r proc~calctdfanddos calcTDFandDOS proc~calctdfanddos->proc~get_hh_r proc~wham_get_eig_deleig wham_get_eig_deleig proc~calctdfanddos->proc~wham_get_eig_deleig proc~wham_get_eig_uu_hh_jjlist wham_get_eig_UU_HH_JJlist proc~wham_get_eig_uu_hh_jjlist->proc~get_hh_r proc~berry_main berry_main proc~berry_main->proc~get_hh_r proc~berry_get_sc_klist berry_get_sc_klist proc~berry_main->proc~berry_get_sc_klist proc~berry_get_shc_klist berry_get_shc_klist proc~berry_main->proc~berry_get_shc_klist proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~berry_main->proc~berry_get_imfgh_klist proc~berry_get_kubo_k berry_get_kubo_k proc~berry_main->proc~berry_get_kubo_k proc~berry_get_imf_klist berry_get_imf_klist proc~berry_main->proc~berry_get_imf_klist proc~k_slice k_slice proc~k_slice->proc~get_hh_r proc~k_slice->proc~wham_get_eig_deleig 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~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~get_hh_r proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_main->proc~gyrotropic_get_k_list proc~wham_get_eig_deleig->proc~get_hh_r proc~wham_get_eig_uu_hh_aa_sc_tb_conv wham_get_eig_UU_HH_AA_sc_TB_conv proc~wham_get_eig_uu_hh_aa_sc_tb_conv->proc~get_hh_r proc~wham_get_eig_uu_hh_aa_sc wham_get_eig_UU_HH_AA_sc proc~wham_get_eig_uu_hh_aa_sc->proc~get_hh_r proc~dos_main dos_main proc~dos_main->proc~get_hh_r proc~dos_main->proc~wham_get_eig_deleig proc~k_path k_path proc~k_path->proc~get_hh_r 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~geninterp_main geninterp_main proc~geninterp_main->proc~get_hh_r proc~spin_get_moment spin_get_moment proc~spin_get_moment->proc~get_hh_r proc~boltzwann_main boltzwann_main proc~boltzwann_main->proc~calctdfanddos proc~berry_get_sc_klist->proc~wham_get_eig_deleig proc~berry_get_sc_klist->proc~wham_get_eig_uu_hh_aa_sc_tb_conv proc~berry_get_sc_klist->proc~wham_get_eig_uu_hh_aa_sc proc~berry_get_shc_klist->proc~wham_get_eig_deleig proc~gyrotropic_get_k_list->proc~wham_get_eig_deleig proc~gyrotropic_get_k_list->proc~berry_get_imfgh_klist proc~gyrotropic_get_k_list->proc~berry_get_imf_klist proc~berry_get_imfgh_klist->proc~wham_get_eig_uu_hh_jjlist proc~berry_get_kubo_k->proc~wham_get_eig_deleig proc~berry_get_imf_klist->proc~berry_get_imfgh_klist

Contents

Source Code


Source Code

  subroutine fourier_q_to_R(op_q, op_R)
    !==========================================================
    !
    !! Fourier transforms Wannier-gauge representation
    !! of a given operator O from q-space to R-space:
    !!
    !! O_ij(q) --> O_ij(R) = (1/N_kpts) sum_q e^{-iqR} O_ij(q)
    !
    !==========================================================

    use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
    use w90_parameters, only: num_kpts, kpt_latt
    use w90_postw90_common, only: nrpts, irvec

    implicit none

    ! Arguments
    !
    complex(kind=dp), dimension(:, :, :), intent(in)  :: op_q
    !! Operator in q-space
    complex(kind=dp), dimension(:, :, :), intent(out) :: op_R
    !! Operator in R-space

    integer          :: ir, ik
    real(kind=dp)    :: rdotq
    complex(kind=dp) :: phase_fac

    op_R = cmplx_0
    do ir = 1, nrpts
      do ik = 1, num_kpts
        rdotq = twopi*dot_product(kpt_latt(:, ik), irvec(:, ir))
        phase_fac = exp(-cmplx_i*rdotq)
        op_R(:, :, ir) = op_R(:, :, ir) + phase_fac*op_q(:, :, ik)
      enddo
    enddo
    op_R = op_R/real(num_kpts, dp)

  end subroutine fourier_q_to_R