spin_get_nk Subroutine

public subroutine spin_get_nk(kpt, spn_nk)

Uses

  • proc~~spin_get_nk~~UsesGraph proc~spin_get_nk spin_get_nk module~w90_postw90_common w90_postw90_common proc~spin_get_nk->module~w90_postw90_common module~w90_utility w90_utility proc~spin_get_nk->module~w90_utility module~w90_constants w90_constants proc~spin_get_nk->module~w90_constants module~w90_get_oper w90_get_oper proc~spin_get_nk->module~w90_get_oper module~w90_io w90_io proc~spin_get_nk->module~w90_io module~w90_parameters w90_parameters proc~spin_get_nk->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 (m=1,...,num_wann) where S.n = n_x.S_x + n_y.S_y + n_z.Z_z

S_i are the Pauli matrices and n=(n_x,n_y,n_z) is the unit vector along the chosen spin quantization axis

Arguments

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

Calls

proc~~spin_get_nk~~CallsGraph proc~spin_get_nk spin_get_nk proc~utility_diagonalize utility_diagonalize proc~spin_get_nk->proc~utility_diagonalize proc~pw90common_fourier_r_to_k pw90common_fourier_R_to_k proc~spin_get_nk->proc~pw90common_fourier_r_to_k proc~utility_rotate_diag utility_rotate_diag proc~spin_get_nk->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

Called by

proc~~spin_get_nk~~CalledByGraph proc~spin_get_nk spin_get_nk proc~k_path k_path proc~k_path->proc~spin_get_nk proc~k_slice k_slice proc~k_slice->proc~spin_get_nk proc~berry_get_kubo_k berry_get_kubo_k proc~berry_get_kubo_k->proc~spin_get_nk proc~berry_main berry_main proc~berry_main->proc~berry_get_kubo_k

Contents

Source Code


Source Code

  subroutine spin_get_nk(kpt, spn_nk)
    !=============================================================!
    !                                                             !
    !! Computes <psi_{mk}^(H)|S.n|psi_{mk}^(H)> (m=1,...,num_wann)
    !! where S.n = n_x.S_x + n_y.S_y + n_z.Z_z
    !!
    !! S_i are the Pauli matrices and n=(n_x,n_y,n_z) is the unit
    !! vector along the chosen spin quantization axis
    !                                                             !
    !============================================================ !

    use w90_constants, only: dp, pi, 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, spin_axis_polar, &
      spin_axis_azimuth
    use w90_postw90_common, only: pw90common_fourier_R_to_k
    use w90_get_oper, only: HH_R, SS_R

    ! Arguments
    !
    real(kind=dp), intent(in)  :: kpt(3)
    real(kind=dp), intent(out) :: spn_nk(num_wann)

    ! Physics
    !
    complex(kind=dp), allocatable :: HH(:, :)
    complex(kind=dp), allocatable :: UU(:, :)
    complex(kind=dp), allocatable :: SS(:, :, :), SS_n(:, :)

    ! Misc/Dummy
    !
    integer          :: is
    real(kind=dp)    :: eig(num_wann), alpha(3), conv

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

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

    do is = 1, 3
      call pw90common_fourier_R_to_k(kpt, SS_R(:, :, :, is), SS(:, :, is), 0)
    enddo

    ! Unit vector along the magnetization direction
    !
    conv = 180.0_dp/pi
    alpha(1) = sin(spin_axis_polar/conv)*cos(spin_axis_azimuth/conv)
    alpha(2) = sin(spin_axis_polar/conv)*sin(spin_axis_azimuth/conv)
    alpha(3) = cos(spin_axis_polar/conv)

    ! Vector of spin matrices projected along the quantization axis
    !
    SS_n(:, :) = alpha(1)*SS(:, :, 1) + alpha(2)*SS(:, :, 2) + alpha(3)*SS(:, :, 3)

    spn_nk(:) = real(utility_rotate_diag(SS_n, UU, num_wann), dp)

  end subroutine spin_get_nk