Given a k point, this function returns eigenvalues E and derivatives of the eigenvalues dE/dk_a, using wham_get_deleig_a
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in), | dimension(3) | :: | kpt | the three coordinates of the k point vector (in relative coordinates) |
|
real(kind=dp), | intent(out) | :: | eig(num_wann) | the calculated eigenvalues at kpt |
||
real(kind=dp), | intent(out) | :: | del_eig(num_wann,3) | the calculated derivatives of the eigenvalues at kpt [first component: band; second component: 1,2,3 for the derivatives along the three k directions] |
||
complex(kind=dp), | intent(out), | dimension(:, :) | :: | HH | the Hamiltonian matrix at kpt |
|
complex(kind=dp), | intent(out), | dimension(:, :, :) | :: | delHH | the delHH matrix (derivative of H) at kpt |
|
complex(kind=dp), | intent(out), | dimension(:, :) | :: | UU | the rotation matrix that gives the eigenvectors of HH |
subroutine wham_get_eig_deleig(kpt, eig, del_eig, HH, delHH, UU)
!! Given a k point, this function returns eigenvalues E and
!! derivatives of the eigenvalues dE/dk_a, using wham_get_deleig_a
!
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
use w90_utility, only: utility_diagonalize
real(kind=dp), dimension(3), intent(in) :: kpt
!! the three coordinates of the k point vector (in relative coordinates)
real(kind=dp), intent(out) :: eig(num_wann)
!! the calculated eigenvalues at kpt
real(kind=dp), intent(out) :: del_eig(num_wann, 3)
!! the calculated derivatives of the eigenvalues at kpt [first component: band; second component: 1,2,3
!! for the derivatives along the three k directions]
complex(kind=dp), dimension(:, :), intent(out) :: HH
!! the Hamiltonian matrix at kpt
complex(kind=dp), dimension(:, :, :), intent(out) :: delHH
!! the delHH matrix (derivative of H) at kpt
complex(kind=dp), dimension(:, :), intent(out) :: UU
!! the rotation matrix that gives the eigenvectors of HH
! I call it to be sure that it has been called already once,
! and that HH_R contains the actual matrix.
! Further calls should return very fast.
call get_HH_R
call pw90common_fourier_R_to_k(kpt, HH_R, HH, 0)
call utility_diagonalize(HH, num_wann, eig, UU)
call pw90common_fourier_R_to_k(kpt, HH_R, delHH(:, :, 1), 1)
call pw90common_fourier_R_to_k(kpt, HH_R, delHH(:, :, 2), 2)
call pw90common_fourier_R_to_k(kpt, HH_R, delHH(:, :, 3), 3)
call wham_get_deleig_a(del_eig(:, 1), eig, delHH(:, :, 1), UU)
call wham_get_deleig_a(del_eig(:, 2), eig, delHH(:, :, 2), UU)
call wham_get_deleig_a(del_eig(:, 3), eig, delHH(:, :, 3), UU)
end subroutine wham_get_eig_deleig