spin_get_moment_k Subroutine

private subroutine spin_get_moment_k(kpt, ef, spn_k)

Uses

  • proc~~spin_get_moment_k~~UsesGraph proc~spin_get_moment_k spin_get_moment_k module~w90_postw90_common w90_postw90_common proc~spin_get_moment_k->module~w90_postw90_common module~w90_utility w90_utility proc~spin_get_moment_k->module~w90_utility module~w90_constants w90_constants proc~spin_get_moment_k->module~w90_constants module~w90_get_oper w90_get_oper proc~spin_get_moment_k->module~w90_get_oper module~w90_io w90_io proc~spin_get_moment_k->module~w90_io module~w90_parameters w90_parameters proc~spin_get_moment_k->module~w90_parameters module~w90_postw90_common->module~w90_constants module~w90_comms w90_comms module~w90_postw90_common->module~w90_comms module~w90_utility->module~w90_constants module~w90_get_oper->module~w90_constants module~w90_io->module~w90_constants module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_comms->module~w90_constants module~w90_comms->module~w90_io

Computes the spin magnetic moment by Wannier interpolation at the specified k-point

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: kpt(3)
real(kind=dp), intent(in) :: ef
real(kind=dp), intent(out) :: spn_k(3)

Calls

proc~~spin_get_moment_k~~CallsGraph proc~spin_get_moment_k spin_get_moment_k proc~utility_diagonalize utility_diagonalize proc~spin_get_moment_k->proc~utility_diagonalize proc~pw90common_get_occ pw90common_get_occ proc~spin_get_moment_k->proc~pw90common_get_occ proc~pw90common_fourier_r_to_k pw90common_fourier_R_to_k proc~spin_get_moment_k->proc~pw90common_fourier_r_to_k proc~utility_rotate_diag utility_rotate_diag proc~spin_get_moment_k->proc~utility_rotate_diag proc~io_error io_error proc~utility_diagonalize->proc~io_error proc~utility_zgemm_new utility_zgemm_new proc~utility_rotate_diag->proc~utility_zgemm_new proc~utility_matmul_diag utility_matmul_diag proc~utility_rotate_diag->proc~utility_matmul_diag

Contents

Source Code


Source Code

  subroutine spin_get_moment_k(kpt, ef, spn_k)
    !! Computes the spin magnetic moment by Wannier interpolation
    !! at the specified k-point
    use w90_constants, only: dp, cmplx_0, cmplx_i
    use w90_io, only: io_error
    use w90_utility, only: utility_diagonalize, utility_rotate_diag
    use w90_parameters, only: num_wann
    use w90_postw90_common, only: pw90common_fourier_R_to_k, pw90common_get_occ
    use w90_get_oper, only: HH_R, SS_R
    ! Arguments
    !
    real(kind=dp), intent(in)  :: kpt(3)
    real(kind=dp), intent(in)  :: ef
    real(kind=dp), intent(out) :: spn_k(3)

    ! Physics
    !
    complex(kind=dp), allocatable :: HH(:, :)
    complex(kind=dp), allocatable :: SS(:, :, :)
    complex(kind=dp), allocatable :: UU(:, :)
    real(kind=dp)                 :: spn_nk(num_wann, 3)

    ! Misc/Dummy
    !
    integer          :: i, is
    real(kind=dp)    :: eig(num_wann), occ(num_wann)

    allocate (HH(num_wann, num_wann))
    allocate (UU(num_wann, num_wann))
    allocate (SS(num_wann, num_wann, 3))

    call pw90common_fourier_R_to_k(kpt, HH_R, HH, 0)
    call utility_diagonalize(HH, num_wann, eig, UU)
    call pw90common_get_occ(eig, occ, ef)

    spn_k(1:3) = 0.0_dp
    do is = 1, 3
      call pw90common_fourier_R_to_k(kpt, SS_R(:, :, :, is), SS(:, :, is), 0)
      spn_nk(:, is) = aimag(cmplx_i*utility_rotate_diag(SS(:, :, is), UU, num_wann))
      do i = 1, num_wann
        spn_k(is) = spn_k(is) + occ(i)*spn_nk(i, is)
      end do
    enddo

  end subroutine spin_get_moment_k