gyrotropic_get_k_list Subroutine

private subroutine gyrotropic_get_k_list(kpt, kweight, gyro_K_spn, gyro_K_orb, gyro_D, gyro_Dw, gyro_C, gyro_DOS, gyro_NOA_orb, gyro_NOA_spn, eval_K, eval_D, eval_Dw, eval_NOA, eval_spn, eval_C, eval_dos)

Uses

  • proc~~gyrotropic_get_k_list~~UsesGraph proc~gyrotropic_get_k_list gyrotropic_get_k_list module~w90_postw90_common w90_postw90_common proc~gyrotropic_get_k_list->module~w90_postw90_common module~w90_utility w90_utility proc~gyrotropic_get_k_list->module~w90_utility module~w90_constants w90_constants proc~gyrotropic_get_k_list->module~w90_constants module~w90_get_oper w90_get_oper proc~gyrotropic_get_k_list->module~w90_get_oper module~w90_spin w90_spin proc~gyrotropic_get_k_list->module~w90_spin module~w90_wan_ham w90_wan_ham proc~gyrotropic_get_k_list->module~w90_wan_ham module~w90_io w90_io proc~gyrotropic_get_k_list->module~w90_io module~w90_parameters w90_parameters proc~gyrotropic_get_k_list->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_spin->module~w90_constants module~w90_wan_ham->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

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in) :: kpt(3)
real(kind=dp), intent(in) :: kweight
real(kind=dp), intent(inout), dimension(:, :, :):: gyro_K_spn
real(kind=dp), intent(inout), dimension(:, :, :):: gyro_K_orb
real(kind=dp), intent(inout), dimension(:, :, :):: gyro_D
real(kind=dp), intent(inout), dimension(:, :, :, :):: gyro_Dw
real(kind=dp), intent(inout), dimension(:, :, :):: gyro_C
real(kind=dp), intent(inout), dimension(:):: gyro_DOS
real(kind=dp), intent(inout), dimension(:, :, :, :):: gyro_NOA_orb
real(kind=dp), intent(inout), dimension(:, :, :, :):: gyro_NOA_spn
logical, intent(in) :: eval_K
logical, intent(in) :: eval_D
logical, intent(in) :: eval_Dw
logical, intent(in) :: eval_NOA
logical, intent(in) :: eval_spn
logical, intent(in) :: eval_C
logical, intent(in) :: eval_dos

Calls

proc~~gyrotropic_get_k_list~~CallsGraph proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~pw90common_fourier_r_to_k_vec pw90common_fourier_R_to_k_vec proc~gyrotropic_get_k_list->proc~pw90common_fourier_r_to_k_vec proc~spin_get_s spin_get_S proc~gyrotropic_get_k_list->proc~spin_get_s proc~gyrotropic_get_noa_k gyrotropic_get_NOA_k proc~gyrotropic_get_k_list->proc~gyrotropic_get_noa_k proc~wham_get_eig_deleig wham_get_eig_deleig proc~gyrotropic_get_k_list->proc~wham_get_eig_deleig proc~wham_get_d_h wham_get_D_h proc~gyrotropic_get_k_list->proc~wham_get_d_h proc~utility_w0gauss utility_w0gauss proc~gyrotropic_get_k_list->proc~utility_w0gauss proc~utility_rotate utility_rotate proc~gyrotropic_get_k_list->proc~utility_rotate proc~berry_get_imf_klist berry_get_imf_klist proc~gyrotropic_get_k_list->proc~berry_get_imf_klist proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~gyrotropic_get_k_list->proc~berry_get_imfgh_klist proc~pw90common_fourier_r_to_k pw90common_fourier_R_to_k proc~spin_get_s->proc~pw90common_fourier_r_to_k proc~utility_rotate_diag utility_rotate_diag proc~spin_get_s->proc~utility_rotate_diag proc~utility_diagonalize utility_diagonalize proc~spin_get_s->proc~utility_diagonalize proc~gyrotropic_get_noa_k->proc~utility_rotate proc~gyrotropic_get_noa_bnl_spin gyrotropic_get_NOA_Bnl_spin proc~gyrotropic_get_noa_k->proc~gyrotropic_get_noa_bnl_spin proc~gyrotropic_get_noa_bnl_orb gyrotropic_get_NOA_Bnl_orb proc~gyrotropic_get_noa_k->proc~gyrotropic_get_noa_bnl_orb proc~get_hh_r get_HH_R proc~wham_get_eig_deleig->proc~get_hh_r proc~wham_get_eig_deleig->proc~pw90common_fourier_r_to_k proc~wham_get_eig_deleig->proc~utility_diagonalize proc~wham_get_d_h->proc~utility_rotate proc~berry_get_imf_klist->proc~berry_get_imfgh_klist proc~berry_get_imfgh_klist->proc~pw90common_fourier_r_to_k_vec proc~wham_get_eig_uu_hh_jjlist wham_get_eig_UU_HH_JJlist proc~berry_get_imfgh_klist->proc~wham_get_eig_uu_hh_jjlist proc~utility_im_tr_prod utility_im_tr_prod proc~berry_get_imfgh_klist->proc~utility_im_tr_prod proc~wham_get_occ_mat_list wham_get_occ_mat_list proc~berry_get_imfgh_klist->proc~wham_get_occ_mat_list proc~utility_re_tr_prod utility_re_tr_prod proc~berry_get_imfgh_klist->proc~utility_re_tr_prod proc~wham_get_eig_uu_hh_jjlist->proc~get_hh_r proc~wham_get_eig_uu_hh_jjlist->proc~utility_diagonalize proc~fourier_q_to_r fourier_q_to_R proc~get_hh_r->proc~fourier_q_to_r proc~io_error io_error proc~get_hh_r->proc~io_error proc~get_win_min get_win_min proc~get_hh_r->proc~get_win_min proc~io_file_unit io_file_unit proc~get_hh_r->proc~io_file_unit proc~utility_matmul_diag utility_matmul_diag proc~utility_rotate_diag->proc~utility_matmul_diag proc~utility_zgemm_new utility_zgemm_new proc~utility_rotate_diag->proc~utility_zgemm_new proc~wham_get_occ_mat_list->proc~io_error proc~utility_diagonalize->proc~io_error

Called by

proc~~gyrotropic_get_k_list~~CalledByGraph proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~gyrotropic_get_k_list

Contents

Source Code


Source Code

  subroutine gyrotropic_get_k_list(kpt, kweight, &
                                   gyro_K_spn, gyro_K_orb, gyro_D, gyro_Dw, gyro_C, &
                                   gyro_DOS, gyro_NOA_orb, gyro_NOA_spn, &
                                   eval_K, eval_D, eval_Dw, eval_NOA, eval_spn, eval_C, eval_dos)
    !======================================================================!
    !                                                                      !
    ! Contribution from point k to the GME tensor, Eq.(9) of ZMS16,        !
    ! evaluated in the clean limit of omega.tau >> 1 where  it is real.    !
    ! The following two quantities are calculated (sigma = Pauli matrix):  !
    !                                                                      !
    ! gyro_K_spn_k = delta(E_kn-E_f).(d E_{kn}/d k_i).sigma_{kn,j}         !
    ! [units of length]                                                    !
    !                                                                      !
    ! gyro_K_orb_k = delta(E_kn-E_f).(d E_{kn}/d k_i).(2.hbar/e).m^orb_{kn,j} !
    ! [units of (length^3)*energy]                                         !
    !                                                                      !
    ! gyro_D_k = delta(E_kn-E_f).(d E_{kn}/d k_i).Omega_{kn,j}             !
    ! [units of length^3]                                                  !
    !                                                                      !
    ! gyro_Dw_k = delta(E_kn-E_f).(d E_{kn}/d k_i).tildeOmega_{kn,j}       !
    ! [units of length^3]                                                  !
    !                                                                      !
    ! gyro_C_k = delta(E_kn-E_f).(d E_{kn}/d k_i).(d E_{kn}/d k_j)         !
    ! [units of energy*length^3]                                           !
    !                                                                      !
    ! gyro_DOS_k = delta(E_kn-E_f)                                          !
    ! [units of 1/Energy]                                                  !
    !                                                                      !
    ! gme_NOA_orb_k = ?????                                           !
    !                                                                      !
    ! gme_NOA_spn_k = ??????                                           !
    !                                                                      !
    !======================================================================!

    use w90_constants, only: dp, cmplx_0, cmplx_i
    use w90_utility, only: utility_rotate, utility_rotate_diag, utility_w0gauss
    use w90_parameters, only: num_wann, fermi_energy_list, &
      gyrotropic_smr_index, nfermi, gyrotropic_nfreq, &
      gyrotropic_degen_thresh, gyrotropic_smr_max_arg, &
      gyrotropic_band_list, gyrotropic_num_bands, &
      gyrotropic_smr_fixed_en_width
    use w90_postw90_common, only: pw90common_get_occ, &
      pw90common_fourier_R_to_k_vec
    use w90_wan_ham, only: wham_get_eig_deleig, wham_get_D_h

    use w90_get_oper, only: HH_R, SS_R, AA_R
    use w90_spin, only: spin_get_S
    use w90_io, only: stdout

    ! Arguments
    !
    real(kind=dp), intent(in)                      :: kpt(3), kweight
    real(kind=dp), dimension(:, :, :), intent(inout)   :: gyro_K_spn, &
                                                          gyro_K_orb, &
                                                          gyro_D, &
                                                          gyro_C
    real(kind=dp), dimension(:, :, :, :), intent(inout) :: gyro_Dw, gyro_NOA_spn, gyro_NOA_orb
    real(kind=dp), dimension(:), intent(inout) :: gyro_DOS

    logical, intent(in) :: eval_K, eval_D, eval_Dw, &
                           eval_C, eval_NOA, eval_spn, eval_dos

    complex(kind=dp), allocatable :: UU(:, :)
    complex(kind=dp), allocatable :: HH(:, :)
    complex(kind=dp), allocatable :: delHH(:, :, :)
    complex(kind=dp), allocatable :: SS(:, :, :)
    complex(kind=dp), allocatable :: AA(:, :, :)
    complex(kind=dp), allocatable :: D_h(:, :, :)

    real(kind=dp), allocatable :: curv_w_nk(:, :, :)

    integer          :: i, j, n, n1, m1, m, ifermi
    real(kind=dp)    :: delta, occ(num_wann), &
                        eig(num_wann), del_eig(num_wann, 3), &
                        S(num_wann, 3), eta_smr, arg, &
                        orb_nk(3), curv_nk(3), &
                        imf_k(3, 3, 1), img_k(3, 3, 1), imh_k(3, 3, 1)
    logical          :: got_spin, got_orb_n

    allocate (UU(num_wann, num_wann))
    allocate (HH(num_wann, num_wann))
    allocate (delHH(num_wann, num_wann, 3))
    allocate (D_h(num_wann, num_wann, 3))

    if (eval_spn) allocate (SS(num_wann, num_wann, 3))

    call wham_get_eig_deleig(kpt, eig, del_eig, HH, delHH, UU)

    if (eval_Dw .or. eval_NOA) then
      allocate (AA(num_wann, num_wann, 3))
      call wham_get_D_h(delHH, UU, eig, D_h)
      call pw90common_fourier_R_to_k_vec(kpt, AA_R, OO_true=AA)
      do i = 1, 3
        AA(:, :, i) = utility_rotate(AA(:, :, i), UU, num_wann)
      enddo
      AA = AA + cmplx_i*D_h ! Eq.(25) WYSV06
    endif

    if (eval_Dw) allocate (curv_w_nk(num_wann, gyrotropic_nfreq, 3))

    eta_smr = gyrotropic_smr_fixed_en_width

    got_spin = .false.

    do n1 = 1, gyrotropic_num_bands
      n = gyrotropic_band_list(n1)
      !
      ! ***ADJUSTABLE PARAMETER***
      ! avoid degeneracies
      !---------------------------------------------------
      if (n > 1) then
        if (eig(n) - eig(n - 1) <= gyrotropic_degen_thresh) cycle
      endif
      if (n < num_wann) then
        if (eig(n + 1) - eig(n) <= gyrotropic_degen_thresh) cycle
      endif
      !---------------------------------------------------
      got_orb_n = .false.
      do ifermi = 1, nfermi
        arg = (eig(n) - fermi_energy_list(ifermi))/eta_smr
        !
        ! To save time: far from the Fermi surface, negligible contribution
        !
        !-------------------------
        if (abs(arg) > gyrotropic_smr_max_arg) cycle
        !-------------------------
        !
        ! Spin is computed for all bands simultaneously
        !
        if (eval_spn .and. .not. got_spin) then
          call spin_get_S(kpt, S)
          got_spin = .true. ! Do it for only one value of ifermi and n
        endif
        ! Orbital quantities are computed for each band separately
        if (.not. got_orb_n) then
          if (eval_K) then
            ! Fake occupations: band n occupied, others empty
            occ = 0.0_dp
            occ(n) = 1.0_dp
            call berry_get_imfgh_klist(kpt, imf_k, img_k, imh_k, occ)
            do i = 1, 3
              orb_nk(i) = sum(imh_k(:, i, 1)) - sum(img_k(:, i, 1))
              curv_nk(i) = sum(imf_k(:, i, 1))
            enddo
          else if (eval_D) then
            occ = 0.0_dp
            occ(n) = 1.0_dp
            call berry_get_imf_klist(kpt, imf_k, occ)
            do i = 1, 3
              curv_nk(i) = sum(imf_k(:, i, 1))
            enddo
            got_orb_n = .true. ! Do it for only one value of ifermi
          endif

          if (eval_Dw) call gyrotropic_get_curv_w_k(eig, AA, curv_w_nk)

          got_orb_n = .true. ! Do it for only one value of ifermi
        endif
        !
        delta = utility_w0gauss(arg, gyrotropic_smr_index)/eta_smr*kweight ! Broadened delta(E_nk-E_f)
        !
        ! Loop over Cartesian tensor components
        !
        do j = 1, 3
          if (eval_K .and. eval_spn) gyro_K_spn(:, j, ifermi) = &
            gyro_K_spn(:, j, ifermi) + del_eig(n, :)*S(n, j)*delta
          if (eval_K) gyro_K_orb(:, j, ifermi) = &
            gyro_K_orb(:, j, ifermi) + del_eig(n, :)*orb_nk(j)*delta
          if (eval_D) gyro_D(:, j, ifermi) = &
            gyro_D(:, j, ifermi) + del_eig(n, :)*curv_nk(j)*delta
          if (eval_Dw) then
            do i = 1, 3
              gyro_Dw(i, j, ifermi, :) = &
                gyro_Dw(i, j, ifermi, :) + del_eig(n, i)*delta*curv_w_nk(n, :, j)
            enddo
          endif
          if (eval_C) gyro_C(:, j, ifermi) = &
            gyro_C(:, j, ifermi) + del_eig(n, :)*del_eig(n, j)*delta
        enddo !j
        if (eval_dos) gyro_DOS(ifermi) = gyro_DOS(ifermi) + delta

      enddo !ifermi
    enddo !n

    if (eval_NOA) then
      if (eval_spn) then
        call gyrotropic_get_NOA_k(kpt, kweight, eig, del_eig, AA, UU, gyro_NOA_orb, gyro_NOA_spn)
      else
        call gyrotropic_get_NOA_k(kpt, kweight, eig, del_eig, AA, UU, gyro_NOA_orb)
      endif
    endif

  end subroutine gyrotropic_get_k_list