wham_get_D_h Subroutine

public subroutine wham_get_D_h(delHH, UU, eig, D_h)

Uses

  • proc~~wham_get_d_h~~UsesGraph proc~wham_get_d_h wham_get_D_h module~w90_utility w90_utility proc~wham_get_d_h->module~w90_utility module~w90_parameters w90_parameters proc~wham_get_d_h->module~w90_parameters module~w90_constants w90_constants proc~wham_get_d_h->module~w90_constants module~w90_utility->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

Compute D^H_a=UU^dag.del_a UU (a=x,y,z) using Eq.(24) of WYSV06

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :, :):: delHH
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
complex(kind=dp), intent(out), dimension(:, :, :):: D_h

Calls

proc~~wham_get_d_h~~CallsGraph proc~wham_get_d_h wham_get_D_h proc~utility_rotate utility_rotate proc~wham_get_d_h->proc~utility_rotate

Called by

proc~~wham_get_d_h~~CalledByGraph proc~wham_get_d_h wham_get_D_h proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->proc~wham_get_d_h proc~berry_get_shc_klist berry_get_shc_klist proc~berry_get_shc_klist->proc~wham_get_d_h proc~berry_get_kubo_k berry_get_kubo_k proc~berry_get_kubo_k->proc~wham_get_d_h proc~k_slice k_slice proc~k_slice->proc~berry_get_shc_klist proc~k_path k_path proc~k_path->proc~berry_get_shc_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_shc_klist proc~berry_main->proc~berry_get_kubo_k

Contents

Source Code


Source Code

  subroutine wham_get_D_h(delHH, UU, eig, D_h)
    !=========================================!
    !                                         !
    !! Compute D^H_a=UU^dag.del_a UU (a=x,y,z)
    !! using Eq.(24) of WYSV06
    !                                         !
    !=========================================!

    ! TO DO: Implement version where energy denominators only connect
    !        occupied and empty states. In this case probably do not need
    !        to worry about avoiding small energy denominators

    use w90_constants, only: dp, cmplx_0
    use w90_parameters, only: num_wann
    use w90_utility, only: utility_rotate

    ! Arguments
    !
    complex(kind=dp), dimension(:, :, :), intent(in)  :: delHH
    complex(kind=dp), dimension(:, :), intent(in)    :: UU
    real(kind=dp), dimension(:), intent(in)    :: eig
    complex(kind=dp), dimension(:, :, :), intent(out) :: D_h

    complex(kind=dp), allocatable :: delHH_bar_i(:, :)
    integer                       :: n, m, i

    allocate (delHH_bar_i(num_wann, num_wann))
    D_h = cmplx_0
    do i = 1, 3
      delHH_bar_i(:, :) = utility_rotate(delHH(:, :, i), UU, num_wann)
      do m = 1, num_wann
        do n = 1, num_wann
          if (n == m .or. abs(eig(m) - eig(n)) < 1.0e-7_dp) cycle
          D_h(n, m, i) = delHH_bar_i(n, m)/(eig(m) - eig(n))
        end do
      end do
    enddo

  end subroutine wham_get_D_h