tran_parity_enforce Subroutine

private subroutine tran_parity_enforce(signatures)

Uses

  • proc~~tran_parity_enforce~~UsesGraph proc~tran_parity_enforce tran_parity_enforce module~w90_constants w90_constants proc~tran_parity_enforce->module~w90_constants module~w90_parameters w90_parameters proc~tran_parity_enforce->module~w90_parameters module~w90_io w90_io proc~tran_parity_enforce->module~w90_io module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_io->module~w90_constants

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(inout), dimension(:, :):: signatures

Called by

proc~~tran_parity_enforce~~CalledByGraph proc~tran_parity_enforce tran_parity_enforce proc~tran_main tran_main proc~tran_main->proc~tran_parity_enforce program~wannier wannier program~wannier->proc~tran_main proc~wannier_run wannier_run proc~wannier_run->proc~tran_main

Contents

Source Code


Source Code

  subroutine tran_parity_enforce(signatures)
    !==============================================================!
    ! Here, the signatures of the each wannier fucntion (stored in !
    ! signatures) is used to determine its relavite parity         !
    ! with respect to the first unit cell. The parity pattern of   !
    ! first unit cell is then enforced.                            !
    !==============================================================!

    use w90_constants, only: dp
    use w90_io, only: stdout, io_stopwatch
    use w90_parameters, only: tran_num_cell_ll, num_wann, tran_num_ll, &
      timing_level, iprint, tran_easy_fix

    implicit none

    real(dp), intent(inout), dimension(:, :)               :: signatures

    integer                                             :: i, j, k, wf_idx, num_wann_cell_ll
    real(dp)                                            :: signature_dot_p

    if (timing_level > 1) call io_stopwatch('tran: parity_enforce', 1)

    !
    ! NP: special "easy" fix of the parities by switching the sign
    ! of the Wannier Functions if the first element of the signature
    ! is found negative. Then updating the signature and the Hamiltonian
    ! matrix element for the corresponding line and column
    !
    if (tran_easy_fix) then
      do i = 1, num_wann
        if (real(signatures(1, i)) .lt. 0.0_dp) then
          signatures(:, i) = -signatures(:, i)
          do k = 1, num_wann
            hr_one_dim(k, i, 0) = -hr_one_dim(k, i, 0)
            hr_one_dim(i, k, 0) = -hr_one_dim(i, k, 0)
          enddo
        endif
      enddo
    endif

    num_wann_cell_ll = tran_num_ll/tran_num_cell_ll
    if (iprint .eq. 5) write (stdout, '(a101)') 'Unit cell    Sorted WF index    Unsort WF index  &
         &Unsorted WF Equiv       Signature Dot Product'
    !
    ! Loop over unit cell in principal layers
    !
    do i = 2, 4*tran_num_cell_ll
      !
      ! Loop over wannier functions in unit cell
      !
      do j = 1, num_wann_cell_ll
        if (i .le. 2*tran_num_cell_ll) then
          wf_idx = j + (i - 1)*num_wann_cell_ll
        else
          wf_idx = num_wann - 2*tran_num_ll + j + (i - 1 - 2*tran_num_cell_ll)*num_wann_cell_ll
        endif
        signature_dot_p = dot_product(signatures(:, tran_sorted_idx(j)), signatures(:, tran_sorted_idx(wf_idx)))
        if (iprint .eq. 5) then
          write (stdout, '(2x,i4,3(13x,i5),12x,f20.17)') &
            i, wf_idx, tran_sorted_idx(wf_idx), tran_sorted_idx(j), signature_dot_p
        endif
        if (abs(signature_dot_p) .le. 0.8_dp) then
          write (stdout, '(a28,i4,a64,i4,a20)') ' WARNING: Wannier function (', tran_sorted_idx(wf_idx), &
            ') seems to has poor resemblance to equivalent wannier function (', tran_sorted_idx(j), ') in first unit cell'
          if (iprint .lt. 5) write (stdout, *) 'Dot product of signatures: ', signature_dot_p
        endif
        if (signature_dot_p .lt. 0.0_dp) then
          do k = 1, num_wann
            hr_one_dim(k, tran_sorted_idx(wf_idx), 0) = -hr_one_dim(k, tran_sorted_idx(wf_idx), 0)
            hr_one_dim(tran_sorted_idx(wf_idx), k, 0) = -hr_one_dim(tran_sorted_idx(wf_idx), k, 0)
          enddo
        endif
      enddo
    enddo

    if (timing_level > 1) call io_stopwatch('tran: parity_enforce', 2)

    return

  end subroutine tran_parity_enforce