pw90common_fourier_R_to_k_vec_dadb_TB_conv Subroutine

public subroutine pw90common_fourier_R_to_k_vec_dadb_TB_conv(kpt, OO_R, OO_da, OO_dadb)

Uses

  • proc~~pw90common_fourier_r_to_k_vec_dadb_tb_conv~~UsesGraph proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv pw90common_fourier_R_to_k_vec_dadb_TB_conv module~w90_ws_distance w90_ws_distance proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv->module~w90_ws_distance module~w90_utility w90_utility proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv->module~w90_utility module~w90_parameters w90_parameters proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv->module~w90_parameters module~w90_constants w90_constants proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv->module~w90_constants module~w90_ws_distance->module~w90_parameters module~w90_ws_distance->module~w90_constants module~w90_utility->module~w90_constants module~w90_parameters->module~w90_constants module~w90_io w90_io module~w90_parameters->module~w90_io module~w90_io->module~w90_constants

For : For :

Arguments

Type IntentOptional AttributesName
real(kind=dp) :: kpt(3)
complex(kind=dp), intent(in), dimension(:, :, :, :):: OO_R
complex(kind=dp), intent(out), optional dimension(:, :, :):: OO_da
complex(kind=dp), intent(out), optional dimension(:, :, :, :):: OO_dadb

Called by

proc~~pw90common_fourier_r_to_k_vec_dadb_tb_conv~~CalledByGraph proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv pw90common_fourier_R_to_k_vec_dadb_TB_conv proc~berry_get_sc_klist berry_get_sc_klist proc~berry_get_sc_klist->proc~pw90common_fourier_r_to_k_vec_dadb_tb_conv proc~berry_main berry_main proc~berry_main->proc~berry_get_sc_klist

Contents


Source Code

  subroutine pw90common_fourier_R_to_k_vec_dadb_TB_conv(kpt, OO_R, OO_da, OO_dadb)
    !====================================================================!
    !                                                                    !
    ! modified version of pw90common_fourier_R_to_k_vec_dadb, includes wannier centres in
    ! the exponential inside the sum (so called TB convention)
    !
    !! For $$OO_{ij;dx,dy,dz}$$:
    !! $$O_{ij;dx,dy,dz}(k) = \sum_R e^{+ik.(R+tau_ij)} O_{ij;dx,dy,dz}(R)$$
    !! For $$OO_{ij;dx1,dy1,dz1;dx2,dy2,dz2}$$:
    !! $$O_{ij;dx1,dy1,dz1;dx2,dy2,dz2}(k) = \sum_R e^{+ik.(R+tau_ij)} i.(R+tau_ij)_{dx2,dy2,dz2}
    !!                                       .O_{ij;dx1,dy1,dz1}(R)$$
    ! with tau_ij = tau_j - tau_i, being tau_i=<0i|r|0i> the individual wannier centres
    !                                                                    !
    !====================================================================!

    use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
    use w90_parameters, only: num_kpts, kpt_latt, num_wann, use_ws_distance, &
      wannier_centres, recip_lattice
    use w90_ws_distance, only: irdist_ws, crdist_ws, wdist_ndeg, ws_translate_dist
    use w90_utility, only: utility_cart_to_frac

    implicit none

    ! Arguments
    !
    real(kind=dp)                                     :: kpt(3)
    complex(kind=dp), dimension(:, :, :, :), intent(in)  :: OO_R
    complex(kind=dp), optional, dimension(:, :, :), intent(out)     :: OO_da
    complex(kind=dp), optional, dimension(:, :, :, :), intent(out)   :: OO_dadb

    integer          :: ir, i, j, ideg, a, b
    real(kind=dp)    :: rdotk
    complex(kind=dp) :: phase_fac
    real(kind=dp)    :: local_wannier_centres(3, num_wann), wannier_centres_frac(3, num_wann)
    real(kind=dp)                                                 :: r_sum(3)

    r_sum = 0.d0

    if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
    if (present(OO_da)) OO_da = cmplx_0
    if (present(OO_dadb)) OO_dadb = cmplx_0

    ! calculate wannier centres in cartesian
    local_wannier_centres(:, :) = 0.d0
    do j = 1, num_wann
      do ir = 1, nrpts
        if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0)) then
          local_wannier_centres(1, j) = real(OO_R(j, j, ir, 1))
          local_wannier_centres(2, j) = real(OO_R(j, j, ir, 2))
          local_wannier_centres(3, j) = real(OO_R(j, j, ir, 3))
        endif
      enddo
    enddo
    ! rotate wannier centres from cartesian to fractional coordinates
    wannier_centres_frac(:, :) = 0.d0
    do ir = 1, num_wann
      call utility_cart_to_frac(local_wannier_centres(:, ir), wannier_centres_frac(:, ir), recip_lattice)
    enddo

!    print *, 'wannier_centres_frac'
!    do ir = 1,num_wann
!     print *, wannier_centres_frac(:,ir)
!    enddo
!    stop
!
!    print *, 'crvec'
!    do ir = 1,nrpts
!     print *, crvec(:,ir)
!    enddo
!    stop
!    print *, 'wannier_centres'
!    do ir = 1,num_wann
!     print *, wannier_centres(:,ir)
!    enddo
!    stop

    do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
      if (use_ws_distance) then
        do j = 1, num_wann
        do i = 1, num_wann
          do ideg = 1, wdist_ndeg(i, j, ir)

            rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir) + &
                                                   wannier_centres_frac(:, j) - wannier_centres_frac(:, i), dp))
            phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
            if (present(OO_da)) then
              ! if we are at the origin and at the same band, then the
              ! matrix element is zero in this convention
              if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0) .and. (i .eq. j)) then
                cycle
              else
                OO_da(i, j, 1) = OO_da(i, j, 1) + phase_fac*OO_R(i, j, ir, 1)
                OO_da(i, j, 2) = OO_da(i, j, 2) + phase_fac*OO_R(i, j, ir, 2)
                OO_da(i, j, 3) = OO_da(i, j, 3) + phase_fac*OO_R(i, j, ir, 3)
              endif
            endif
            if (present(OO_dadb)) then
              ! same skip as before
              if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0) .and. (i .eq. j)) then
                cycle
              else
                do a = 1, 3
                  do b = 1, 3
                    OO_dadb(i, j, a, b) = OO_dadb(i, j, a, b) + cmplx_i* &
                                          (crdist_ws(b, ideg, i, j, ir) + local_wannier_centres(b, j) - &
                                           local_wannier_centres(b, i))*phase_fac*OO_R(i, j, ir, a)
                  enddo
                enddo
              endif
            endif

          enddo
        enddo
        enddo
      else
! [lp] Original code, without IJ-dependent shift:
        do j = 1, num_wann
        do i = 1, num_wann
          r_sum(:) = real(irvec(:, ir)) + wannier_centres_frac(:, j) - wannier_centres_frac(:, i)
          rdotk = twopi*dot_product(kpt(:), r_sum(:))
          phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir), dp)
          if (present(OO_da)) then
            ! if we are at the origin and at the same band, then the
            ! matrix element is zero in this convention
            if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0) .and. (i .eq. j)) then
              OO_da(i, j, 1) = OO_da(i, j, 1) + phase_fac*(OO_R(i, j, ir, 1) - local_wannier_centres(1, j))
              OO_da(i, j, 2) = OO_da(i, j, 2) + phase_fac*(OO_R(i, j, ir, 2) - local_wannier_centres(2, j))
              OO_da(i, j, 3) = OO_da(i, j, 3) + phase_fac*(OO_R(i, j, ir, 3) - local_wannier_centres(3, j))
!            print *, 'OO_R(i,j,ir,1)', OO_R(i,j,ir,1)
!            print *, 'local_wannier_centres(1,j)', local_wannier_centres(1,j)
!            print *, 'OO_R(i,j,ir,2)', OO_R(i,j,ir,2)
!            print *, 'local_wannier_centres(2,j)', local_wannier_centres(2,j)
              cycle
            else
              OO_da(i, j, 1) = OO_da(i, j, 1) + phase_fac*OO_R(i, j, ir, 1)
              OO_da(i, j, 2) = OO_da(i, j, 2) + phase_fac*OO_R(i, j, ir, 2)
              OO_da(i, j, 3) = OO_da(i, j, 3) + phase_fac*OO_R(i, j, ir, 3)
            endif
          endif
          if (present(OO_dadb)) then
            ! same skip as before
            if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0) .and. (i .eq. j)) then
              do a = 1, 3
                do b = 1, 3
                  OO_dadb(i, j, a, b) = OO_dadb(i, j, a, b) + cmplx_i*(crvec(b, ir) + local_wannier_centres(b, j) - &
                                                                       local_wannier_centres(b, i))*phase_fac* &
                                        (OO_R(i, j, ir, a) - local_wannier_centres(a, j))
                enddo
              enddo
!           cycle
            else
              do a = 1, 3
                do b = 1, 3
                  OO_dadb(i, j, a, b) = OO_dadb(i, j, a, b) + cmplx_i*(crvec(b, ir) + local_wannier_centres(b, j) - &
                                                                       local_wannier_centres(b, i))*phase_fac*OO_R(i, j, ir, a)
                enddo
              enddo
            endif
          endif
        enddo
        enddo
      endif
    enddo

  end subroutine pw90common_fourier_R_to_k_vec_dadb_TB_conv