Occupation matrix f, and g=1-f for a list of Fermi energies
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=dp), | intent(in), | dimension(:, :) | :: | UU | ||
complex(kind=dp), | intent(out), | dimension(:, :, :) | :: | f_list | ||
complex(kind=dp), | intent(out), | dimension(:, :, :) | :: | g_list | ||
real(kind=dp), | intent(in), | optional | dimension(:) | :: | eig | |
real(kind=dp), | intent(in), | optional | dimension(:) | :: | occ |
subroutine wham_get_occ_mat_list(UU, f_list, g_list, eig, occ)
! subroutine wham_get_occ_mat_list(eig,UU,f_list,g_list)
!================================!
! !
!! Occupation matrix f, and g=1-f
!! for a list of Fermi energies
! Tsirkin: !now optionally either eig or occ parameters may be supplied !
! (Changed consistently the calls from the Berry module) !
!================================!
use w90_constants, only: dp, cmplx_0, cmplx_1
use w90_parameters, only: num_wann, nfermi, fermi_energy_list
use w90_postw90_common, only: pw90common_get_occ
use w90_io, only: io_error
! Arguments
!
complex(kind=dp), dimension(:, :), intent(in) :: UU
complex(kind=dp), dimension(:, :, :), intent(out) :: f_list
complex(kind=dp), dimension(:, :, :), intent(out) :: g_list
real(kind=dp), intent(in), optional, dimension(:) :: eig
real(kind=dp), intent(in), optional, dimension(:) :: occ
integer :: n, m, i, if, nfermi_loc
real(kind=dp), allocatable :: occ_list(:, :)
if (present(occ)) then
nfermi_loc = 1
else
nfermi_loc = nfermi
endif
allocate (occ_list(num_wann, nfermi_loc))
if (present(occ) .and. present(eig)) then
call io_error( &
'occ_list and eig cannot be both arguments in get_occ_mat_list')
elseif (.not. present(occ) .and. .not. present(eig)) then
call io_error( &
'either occ_list or eig must be passed as arguments to get_occ_mat_list')
endif
if (present(occ)) then
occ_list(:, 1) = occ(:)
else
do if = 1, nfermi_loc
call pw90common_get_occ(eig, occ_list(:, if), fermi_energy_list(if))
enddo
endif
f_list = cmplx_0
do if = 1, nfermi_loc
do n = 1, num_wann
do m = 1, num_wann
do i = 1, num_wann
f_list(n, m, if) = f_list(n, m, if) &
+ UU(n, i)*occ_list(i, if)*conjg(UU(m, i))
enddo
g_list(n, m, if) = -f_list(n, m, if)
if (m == n) g_list(n, n, if) = g_list(n, n, if) + cmplx_1
enddo
enddo
enddo
end subroutine wham_get_occ_mat_list