wham_get_D_h_P_value Subroutine

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

Uses

  • proc~~wham_get_d_h_p_value~~UsesGraph proc~wham_get_d_h_p_value wham_get_D_h_P_value module~w90_utility w90_utility proc~wham_get_d_h_p_value->module~w90_utility module~w90_parameters w90_parameters proc~wham_get_d_h_p_value->module~w90_parameters module~w90_constants w90_constants proc~wham_get_d_h_p_value->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_p_value~~CallsGraph proc~wham_get_d_h_p_value wham_get_D_h_P_value proc~utility_rotate utility_rotate proc~wham_get_d_h_p_value->proc~utility_rotate

Called by

proc~~wham_get_d_h_p_value~~CalledByGraph proc~wham_get_d_h_p_value wham_get_D_h_P_value proc~berry_get_sc_klist berry_get_sc_klist proc~berry_get_sc_klist->proc~wham_get_d_h_p_value proc~berry_main berry_main proc~berry_main->proc~berry_get_sc_klist

Contents

Source Code


Source Code

  subroutine wham_get_D_h_P_value(delHH, UU, eig, D_h)
    !=========================================!
    !                                         !
    !! Compute D^H_a=UU^dag.del_a UU (a=x,y,z)
    !! using Eq.(24) of WYSV06
    !  and prescription for energy denominator
    !  from BK81
    !                                         !
    !=========================================!

    ! 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, sc_eta
    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
    real(kind=dp)                 :: deltaE

    allocate (delHH_bar_i(num_wann, num_wann))
    D_h = cmplx_0
    deltaE = 0.d0
    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) cycle
          deltaE = eig(m) - eig(n)
          D_h(n, m, i) = delHH_bar_i(n, m)*(deltaE/(deltaE**(2) + sc_eta**(2)))
        end do
      end do
    enddo

  end subroutine wham_get_D_h_P_value