Find the supercell translation (i.e. the translation by a integer number of supercell vectors, the supercell being defined by the mp_grid) that minimizes the distance between two given Wannier functions, i and j, the first in unit cell 0, the other in unit cell R. I.e., we find the translation to put WF j in the Wigner-Seitz of WF i. We also look for the number of equivalent translation, that happen when w_j,R is on the edge of the WS of w_i,0. The results are stored in global arrays wdist_ndeg, irdist_ws, crdist_ws.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nrpts | |||
integer, | intent(in) | :: | irvec(3,nrpts) | |||
logical, | intent(in), | optional | :: | force_recompute |
subroutine ws_translate_dist(nrpts, irvec, force_recompute)
!! Find the supercell translation (i.e. the translation by a integer number of
!! supercell vectors, the supercell being defined by the mp_grid) that
!! minimizes the distance between two given Wannier functions, i and j,
!! the first in unit cell 0, the other in unit cell R.
!! I.e., we find the translation to put WF j in the Wigner-Seitz of WF i.
!! We also look for the number of equivalent translation, that happen when w_j,R
!! is on the edge of the WS of w_i,0. The results are stored in global
!! arrays wdist_ndeg, irdist_ws, crdist_ws.
use w90_parameters, only: num_wann, wannier_centres, real_lattice, &
recip_lattice, iprint
!translation_centre_frac, automatic_translation,lenconfac
use w90_io, only: stdout, io_error
use w90_utility, only: utility_cart_to_frac, utility_frac_to_cart
implicit none
integer, intent(in) :: nrpts
integer, intent(in) :: irvec(3, nrpts)
logical, optional, intent(in):: force_recompute ! set to true to force recomputing everything
! <<<local variables>>>
integer :: iw, jw, ideg, ir, ierr
integer :: shifts(3, ndegenx)
real(DP) :: irvec_cart(3), tmp(3), tmp_frac(3), R_out(3, ndegenx)
! The subroutine does nothing if called more than once, which may
! not be the best thing if you invoke it while the WFs are moving
if (present(force_recompute)) then
if (force_recompute) then
call clean_ws_translate()
endif
endif
if (done_ws_distance) return
done_ws_distance = .true.
if (ndegenx*num_wann*nrpts <= 0) then
call io_error("unexpected dimensions in ws_translate_dist")
end if
allocate (irdist_ws(3, ndegenx, num_wann, num_wann, nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating irdist_ws in ws_translate_dist')
allocate (crdist_ws(3, ndegenx, num_wann, num_wann, nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating crdist_ws in ws_translate_dist')
allocate (wdist_ndeg(num_wann, num_wann, nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating wcenter_ndeg in ws_translate_dist')
!translation_centre_frac = 0._dp
wdist_ndeg = 0
irdist_ws = 0
crdist_ws = 0
do ir = 1, nrpts
do jw = 1, num_wann
do iw = 1, num_wann
call utility_frac_to_cart(REAL(irvec(:, ir), kind=dp), irvec_cart, real_lattice)
! function JW translated in the Wigner-Seitz around function IW
! and also find its degeneracy, and the integer shifts needed
! to identify it
! Note: the routine outputs R_out, but we don't really need it
! This is kept in case in the future we might want to use it
! R_out contains the actual vector between the two WFs. We
! calculate instead crdist_ws, that is the Bravais lattice vector
! between two supercell lattices, that is the only one we need
! later for interpolation etc.
CALL R_wz_sc(-wannier_centres(:, iw) &
+ (irvec_cart + wannier_centres(:, jw)), (/0._dp, 0._dp, 0._dp/), &
wdist_ndeg(iw, jw, ir), R_out, shifts)
do ideg = 1, wdist_ndeg(iw, jw, ir)
irdist_ws(:, ideg, iw, jw, ir) = irvec(:, ir) + shifts(:, ideg)
tmp_frac = REAL(irdist_ws(:, ideg, iw, jw, ir), kind=dp)
CALL utility_frac_to_cart(tmp_frac, tmp, real_lattice)
crdist_ws(:, ideg, iw, jw, ir) = tmp
enddo
enddo
enddo
enddo
end subroutine ws_translate_dist