wham_get_D_h_a Subroutine

private subroutine wham_get_D_h_a(delHH_a, UU, eig, ef, D_h_a)

Uses

  • proc~~wham_get_d_h_a~~UsesGraph proc~wham_get_d_h_a wham_get_D_h_a module~w90_postw90_common w90_postw90_common proc~wham_get_d_h_a->module~w90_postw90_common module~w90_utility w90_utility proc~wham_get_d_h_a->module~w90_utility module~w90_parameters w90_parameters proc~wham_get_d_h_a->module~w90_parameters module~w90_constants w90_constants proc~wham_get_d_h_a->module~w90_constants module~w90_postw90_common->module~w90_constants module~w90_comms w90_comms module~w90_postw90_common->module~w90_comms 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_comms->module~w90_constants module~w90_comms->module~w90_io module~w90_io->module~w90_constants

Compute D^H_a=UU^dag.del_a UU (a=alpha,beta), using Eq.(24) of WYSV06

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :):: delHH_a
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
real(kind=dp), intent(in) :: ef
complex(kind=dp), intent(out), dimension(:, :):: D_h_a

Calls

proc~~wham_get_d_h_a~~CallsGraph proc~wham_get_d_h_a wham_get_D_h_a proc~pw90common_get_occ pw90common_get_occ proc~wham_get_d_h_a->proc~pw90common_get_occ proc~utility_rotate utility_rotate proc~wham_get_d_h_a->proc~utility_rotate

Contents

Source Code


Source Code

  subroutine wham_get_D_h_a(delHH_a, UU, eig, ef, D_h_a)
    !===============================================!
    !                                               !
    !! Compute D^H_a=UU^dag.del_a UU (a=alpha,beta),
    !! using Eq.(24) of WYSV06
    !                                               !
    !===============================================!

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

    ! Arguments
    !
    complex(kind=dp), dimension(:, :), intent(in)  :: delHH_a
    complex(kind=dp), dimension(:, :), intent(in)  :: UU
    real(kind=dp), dimension(:), intent(in)  :: eig
    real(kind=dp), intent(in)  :: ef
    complex(kind=dp), dimension(:, :), intent(out) :: D_h_a

    complex(kind=dp), allocatable :: delHH_a_bar(:, :)
    real(kind=dp)                 :: occ(num_wann)
    integer                       :: n, m

    call pw90common_get_occ(eig, occ, ef)

    allocate (delHH_a_bar(num_wann, num_wann))
    delHH_a_bar = utility_rotate(delHH_a, UU, num_wann)
    do m = 1, num_wann
      do n = 1, num_wann
        if (occ(n) > 0.999_dp .and. occ(m) < 0.001_dp) then
          D_h_a(n, m) = delHH_a_bar(n, m)/(eig(m) - eig(n))
        else
          D_h_a(n, m) = cmplx_0
        end if
      end do
    end do
    D_h_a = D_h_a - conjg(transpose(D_h_a))

  end subroutine wham_get_D_h_a