Compute D^H_a=UU^dag.del_a UU (a=alpha,beta), using Eq.(24) of WYSV06
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
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