wham_get_eig_UU_HH_JJlist Subroutine

public subroutine wham_get_eig_UU_HH_JJlist(kpt, eig, UU, HH, JJp_list, JJm_list, occ)

Uses

  • proc~~wham_get_eig_uu_hh_jjlist~~UsesGraph proc~wham_get_eig_uu_hh_jjlist wham_get_eig_UU_HH_JJlist module~w90_postw90_common w90_postw90_common proc~wham_get_eig_uu_hh_jjlist->module~w90_postw90_common module~w90_get_oper w90_get_oper proc~wham_get_eig_uu_hh_jjlist->module~w90_get_oper module~w90_utility w90_utility proc~wham_get_eig_uu_hh_jjlist->module~w90_utility module~w90_parameters w90_parameters proc~wham_get_eig_uu_hh_jjlist->module~w90_parameters module~w90_comms w90_comms module~w90_postw90_common->module~w90_comms module~w90_constants w90_constants module~w90_postw90_common->module~w90_constants module~w90_get_oper->module~w90_constants module~w90_utility->module~w90_constants module~w90_io w90_io module~w90_parameters->module~w90_io module~w90_parameters->module~w90_constants module~w90_io->module~w90_constants module~w90_comms->module~w90_io module~w90_comms->module~w90_constants

Wrapper routine used to reduce number of Fourier calls

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt
real(kind=dp), intent(out) :: eig(num_wann)
complex(kind=dp), intent(out), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :):: HH
complex(kind=dp), intent(out), dimension(:, :, :, :):: JJp_list
complex(kind=dp), intent(out), dimension(:, :, :, :):: JJm_list
real(kind=dp), intent(in), optional dimension(:):: occ

Calls

proc~~wham_get_eig_uu_hh_jjlist~~CallsGraph proc~wham_get_eig_uu_hh_jjlist wham_get_eig_UU_HH_JJlist proc~utility_diagonalize utility_diagonalize proc~wham_get_eig_uu_hh_jjlist->proc~utility_diagonalize proc~get_hh_r get_HH_R proc~wham_get_eig_uu_hh_jjlist->proc~get_hh_r proc~io_error io_error proc~utility_diagonalize->proc~io_error proc~get_win_min get_win_min proc~get_hh_r->proc~get_win_min proc~fourier_q_to_r fourier_q_to_R proc~get_hh_r->proc~fourier_q_to_r proc~get_hh_r->proc~io_error proc~io_file_unit io_file_unit proc~get_hh_r->proc~io_file_unit

Called by

proc~~wham_get_eig_uu_hh_jjlist~~CalledByGraph proc~wham_get_eig_uu_hh_jjlist wham_get_eig_UU_HH_JJlist proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~berry_get_imfgh_klist->proc~wham_get_eig_uu_hh_jjlist proc~berry_get_imf_klist berry_get_imf_klist proc~berry_get_imf_klist->proc~berry_get_imfgh_klist proc~k_slice k_slice 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_imfgh_klist proc~k_path->proc~berry_get_imf_klist proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->proc~berry_get_imfgh_klist proc~gyrotropic_get_k_list->proc~berry_get_imf_klist proc~berry_main berry_main proc~berry_main->proc~berry_get_imfgh_klist proc~berry_main->proc~berry_get_imf_klist proc~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~gyrotropic_get_k_list

Contents


Source Code

  subroutine wham_get_eig_UU_HH_JJlist(kpt, eig, UU, HH, JJp_list, JJm_list, occ)
    !========================================================!
    !                                                        !
    !! Wrapper routine used to reduce number of Fourier calls
    !    Added the optional occ parameter                    !
    !========================================================!

    use w90_parameters, only: num_wann
    use w90_get_oper, only: HH_R, get_HH_R
    use w90_postw90_common, only: pw90common_fourier_R_to_k_new
    use w90_utility, only: utility_diagonalize

    real(kind=dp), dimension(3), intent(in)           :: kpt
    real(kind=dp), intent(out)                        :: eig(num_wann)
    complex(kind=dp), dimension(:, :), intent(out)     :: UU
    complex(kind=dp), dimension(:, :), intent(out)     :: HH
    complex(kind=dp), dimension(:, :, :, :), intent(out) :: JJp_list
    complex(kind=dp), dimension(:, :, :, :), intent(out) :: JJm_list
    real(kind=dp), intent(in), optional, dimension(:) :: occ

    integer                       :: i
    complex(kind=dp), allocatable :: delHH(:, :, :)

    call get_HH_R

    allocate (delHH(num_wann, num_wann, 3))
    call pw90common_fourier_R_to_k_new(kpt, HH_R, OO=HH, &
                                       OO_dx=delHH(:, :, 1), &
                                       OO_dy=delHH(:, :, 2), &
                                       OO_dz=delHH(:, :, 3))
    call utility_diagonalize(HH, num_wann, eig, UU)
    do i = 1, 3
      if (present(occ)) then
        call wham_get_JJp_JJm_list(delHH(:, :, i), UU, eig, JJp_list(:, :, :, i), JJm_list(:, :, :, i), occ=occ)
      else
        call wham_get_JJp_JJm_list(delHH(:, :, i), UU, eig, JJp_list(:, :, :, i), JJm_list(:, :, :, i))
      endif
    enddo

  end subroutine wham_get_eig_UU_HH_JJlist